English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Optimizer: TSumOfSquares.h Source File
Uranie / Optimizer  v4.10.0
/* @license-end */
TSumOfSquares.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 
35 #ifndef TSUMOFSQUARES_H
36 #define TSUMOFSQUARES_H
37 
38 // ROOT
39 #include "TH1.h"
40 #include "TF1.h"
41 #include "TPad.h"
42 #include "TRandom3.h"
43 #include "TVirtualFitter.h"
44 #include "TStyle.h"
45 #include "Minuit2/FCNBase.h"
46 #include "TMatrixT.h"
47 
48 #include <vector>
49 #include <iostream>
50 
51 // URANIE
52 #include "Optimizer.h"
53 #include "TCode.h"
54 #include "TDataServer.h"
55 #include "TObjective.h"
56 
57 namespace URANIE
58 {
59 namespace Optimizer
60 {
61 class TSumOfSquares: public ROOT::Minuit2::FCNBase
62 {
63 private:
64  vector<TObjective *> _objectives;
65 
66 public:
67 
69  _objectives(0), _ninput(0), _noutput(0), _tds(NULL), _code(NULL), _sCriteria("")
70  {
71  }
73  TSumOfSquares(URANIE::DataServer::TDataServer *tds,
74  URANIE::Launcher::TCode *code) :
75  _objectives(0), _tds(tds), _code(code), _sCriteria("")
76  {
77  init();
78  }
79 
81  {
82  _dinputValue = NULL;
83  delete[] _dinputValue;
84  _doutputValue = NULL;
85  delete[] _doutputValue;
86 
87  }
88 
89  void init();
90 
91  double operator()(const vector<double> & par) const;
92 
93  double Up() const
94  {
95  return 1.;
96  }
97 
98  void clean();
99 
101  {
102  _objectives.push_back(tobj);
103  }
105  {
106  return (Int_t)(_objectives.size());
107  }
108 
109 private:
110  int _ninput;
111  double * _dinputValue;
112  int _noutput;
113  double * _doutputValue;
114  URANIE::DataServer::TDataServer *_tds;
115  URANIE::Launcher::TCode *_code;
117  TString _sCriteria;
118 
119  ClassDef(URANIE::Optimizer::TSumOfSquares, ID_OPTIMIZER)
120  //Classe de
121 };
122 } // Fin du namespace Optimizer
123 } // Fin du namespace URANIE
124 #endif
125 
Rosenbrock&#39;s function (n=2) with first and second order derivatives.
Definition: TBestEstimate.h:57
TString _sCriteria
List of matrix.
Definition: TSumOfSquares.h:117
int _noutput
Definition: TSumOfSquares.h:112
TList * _listOfDataOutput
Pointer vers un TCode.
Definition: TSumOfSquares.h:116
Int_t getNObjectives()
Definition: TSumOfSquares.h:104
void addObjective(TObjective *tobj)
Definition: TSumOfSquares.h:100
URANIE::Launcher::TCode * _code
Pointer vers un TDS.
Definition: TSumOfSquares.h:115
TSumOfSquares(URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *code)
Constructor with a dataserver.
Definition: TSumOfSquares.h:73
double * _doutputValue
Definition: TSumOfSquares.h:113
int _ninput
Definition: TSumOfSquares.h:110
double Up() const
Definition: TSumOfSquares.h:93
URANIE::DataServer::TDataServer * _tds
Definition: TSumOfSquares.h:114
~TSumOfSquares()
Definition: TSumOfSquares.h:80
vector< TObjective * > _objectives
vector of objectives
Definition: TSumOfSquares.h:64
Interface of class URANIE::Optimize::TObjective.
TSumOfSquares()
Definition: TSumOfSquares.h:68
double operator()(const vector< double > &par) const
Definition: TSumOfSquares.h:61
Description of the class TObjective. This class computes the ojective (from L2 point of view)...
Definition: TObjective.h:61
double * _dinputValue
Definition: TSumOfSquares.h:111