English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Optimizer: TMultiGenSumOfSquares.h Source File
Uranie / Optimizer v4.9.0
/* @license-end */
TMultiGenSumOfSquares.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
38#ifndef TMULTIGENSUMOFSQUARES_H
39#define TMULTIGENSUMOFSQUARES_H
40
41// ROOT
42#include "Math/IFunction.h"
43#include "TH1.h"
44#include "TF1.h"
45#include "TFile.h"
46#include "TPad.h"
47#include "TRandom3.h"
48#include "TStyle.h"
49#include "TSystem.h"
50
51#include <vector>
52#include <iostream>
53
54// URANIE
55#include "Optimizer.h"
56#include "TCode.h"
57#include "TDataServer.h"
58#include "TObjective.h"
59
60namespace URANIE
61{
62namespace Optimizer
63{
64class TMultiGenSumOfSquares: public ROOT::Math::IBaseFunctionMultiDim
65{
66private:
67 vector<TObjective *> _objectives;
68public:
69 //---------------------------------------------
73
74 TMultiGenSumOfSquares(URANIE::DataServer::TDataServer *tds,
75 URANIE::Launcher::TCode *code);
78
79 //---------------------------------------------
83
84 void init();
86 void clean();
90 {
91 Int_t nres = (Int_t)(_objectives.size());
92 return nres;
93 }
94 const char * getCriteriaName()
95 {
96 return _sCriteria.Data();
97 }
98 Double_t getSumOfWeight()
99 {
100 return _dSumOfWeight;
101 }
102 //---------------------------------------------
106
108 {
110 _code);
111 for (Int_t indx = 0; indx < (Int_t)(_objectives.size()); indx++)
112 theClone->addObjective(_objectives[indx]);
113 theClone->init();
114 return theClone;
115 }
116
118 unsigned int NDim() const
119 {
120 return _ninput;
121 }
123
124
125 //---------------------------------------------
129 void setLog()
130 {
131 _blog = kTRUE;
132 }
133 void unsetLog()
134 {
135 _blog = kFALSE;
136 }
137
139
140
141private:
142 //---------------------------------------------
146
147 double DoEval(const double * x) const;
149
150 Bool_t _blog;
152 double * _dinputValue;
155 Double_t _dSumOfWeight;
156 URANIE::DataServer::TDataServer *_tds;
157 URANIE::Launcher::TCode *_code;
159 TString _sCriteria;
160
161 ClassDef(URANIE::Optimizer::TMultiGenSumOfSquares, ID_OPTIMIZER)
162 //Classe de
163
164};
165} // Fin du namespace Optimizer
166} // Fin du namespace URANIE
167
169
170#endif
R__EXTERN URANIE::Optimizer::TMultiGenSumOfSquares * gTMGSOS
Definition TMultiGenSumOfSquares.h:168
Interface of class URANIE::Optimize::TObjective.
NEWMAT::ColumnVector & x
Definition TOptimizerOpt.cxx:70
Definition TMultiGenSumOfSquares.h:65
Double_t _dSumOfWeight
The sum of weights.
Definition TMultiGenSumOfSquares.h:155
Int_t getNObjectives()
Definition TMultiGenSumOfSquares.h:89
double * _dinputValue
Definition TMultiGenSumOfSquares.h:152
int _ninput
Definition TMultiGenSumOfSquares.h:151
double DoEval(const double *x) const
the DoEval method
int _noutput
Definition TMultiGenSumOfSquares.h:153
Bool_t _blog
Boolean for edit the log.
Definition TMultiGenSumOfSquares.h:150
const char * getCriteriaName()
Definition TMultiGenSumOfSquares.h:94
TMultiGenSumOfSquares * Clone() const
The clone method.
Definition TMultiGenSumOfSquares.h:107
void unsetLog()
Definition TMultiGenSumOfSquares.h:133
double * _doutputValue
Definition TMultiGenSumOfSquares.h:154
URANIE::DataServer::TDataServer * _tds
Definition TMultiGenSumOfSquares.h:156
TList * _listOfDataOutput
Pointer vers un TCode.
Definition TMultiGenSumOfSquares.h:158
Double_t getSumOfWeight()
Definition TMultiGenSumOfSquares.h:98
unsigned int NDim() const
the NDIM method
Definition TMultiGenSumOfSquares.h:118
void setLog()
Definition TMultiGenSumOfSquares.h:129
TMultiGenSumOfSquares(URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *code)
Constructor with a dataserver.
TString _sCriteria
List of matrix.
Definition TMultiGenSumOfSquares.h:159
vector< TObjective * > _objectives
vector of objectives
Definition TMultiGenSumOfSquares.h:67
URANIE::Launcher::TCode * _code
Pointer vers un TDS.
Definition TMultiGenSumOfSquares.h:157
Description of the class TObjective. This class computes the ojective (from L2 point of view).
Definition TObjective.h:62
Rosenbrock's function (n=2) with first and second order derivatives.
Definition TBestEstimate.h:58