English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Modeler: TLinearRegression.h Source File
Uranie / Modeler v4.9.0
/* @license-end */
TLinearRegression.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
41#ifndef TLINEARREGRESSION_H
42#define TLINEARREGRESSION_H
43#include "TModeler.h"
44
45#include <iostream>
46
47namespace URANIE
48{
49namespace Modeler
50{
52{
53
54private:
55 TMatrixD computeH(TMatrixD A);
56 TMatrixD computeCov(TMatrixD A);
57public:
58 TVectorD computeHii(TMatrixD A);
59public:
60 Double_t _dr2a;
61 Double_t _dq2;
62 Double_t _dpress;
64
65 // Operations
66public:
67 //---------------------------------------------
71
72 /*
73 \param tds (URANIE::DataServer::TDataServer *) the dataserver
74 \param varexpinput(const char*)[""] The list of input attributes to pass to the function separated by the caracter ":"
75 \param varexpoutput (const char *)[""] The list of output attributes separated by the caracter ":"
76 \param
77 */
78 TLinearRegression(URANIE::DataServer::TDataServer *tds,
79 TString architecture, Option_t *option = "");
81 TLinearRegression(URANIE::DataServer::TDataServer *tds,
82 const char *varexpinput, const char *varexpoutput,
83 Option_t * option = "");
85#ifdef WITH_LIBXML2
86 TLinearRegression (URANIE::DataServer::TDataServer *tds,
87 int id,
88 const char* pmmlfile,
89 const char* LRname,
90 Option_t * option="");
91#endif
95
96 //---------------------------------------------
100
101 Double_t getR2Adjusted() const
102 {
103 return _dr2a;
104 }
106 Double_t getPress() const
107 {
108 return _dpress;
109 }
111 Double_t getQ2() const
112 {
113 return _dq2;
114 }
116 TMatrixD getParametersCovariance() const
117 {
118 return _mParameterCov;
119 }
121 //---------------------------------------------
125 void estimate(Option_t * option = "");
127
128 //---------------------------------------------
132
137 void exportModelCplusplus(std::ofstream * sourcefile) const;
139
143 void exportModelFortran(std::ofstream * sourcefile) const;
144
146 void exportModelPMML(const char* file = "", const char* name = "",
147 Option_t *option = "") const;
148
149
151
155 void exportModelPython(std::ofstream * sourcefile) const
156 {
157 UNUSED(sourcefile);
158 std::cerr
159 << " --- Method TLinearRegression::exportModelPython(ofstream * sourcefile) not yet implemented"
160 << endl;
161 };
163
164 //---------------------------------------------
168 virtual void printLog(Option_t *option = "");
170
171 ClassDef(URANIE::Modeler::TLinearRegression, ID_MODELER)
172 //Classe de
173};
174
175} // Fin du namespace Modeler
176} // Fin du namespace URANIE
177#endif
Interface of the class URANIE::Optimize::TModeler.
Description of the class TLinearRegression.
Definition TLinearRegression.h:52
void estimate(Option_t *option="")
Double_t _dpress
The PRESS quality.
Definition TLinearRegression.h:62
Double_t getQ2() const
Return the Q2 statistic.
Definition TLinearRegression.h:111
Double_t getPress() const
Return the Press statistic.
Definition TLinearRegression.h:106
Double_t _dq2
The Q2 quality.
Definition TLinearRegression.h:61
void exportModelPMML(const char *file="", const char *name="", Option_t *option="") const
Export the model in a PMML file.
void exportModelCplusplus(std::ofstream *sourcefile) const
Export the model in C++ langage in a file.
TMatrixD _mParameterCov
Definition TLinearRegression.h:63
TLinearRegression(URANIE::DataServer::TDataServer *tds, const char *varexpinput, const char *varexpoutput, Option_t *option="")
Default constructor with input and output attributes.
TLinearRegression(URANIE::DataServer::TDataServer *tds, TString architecture, Option_t *option="")
Default constructor with the TDataServer.
TMatrixD computeCov(TMatrixD A)
Double_t getR2Adjusted() const
Return the Adjusted R2 statistic.
Definition TLinearRegression.h:101
TVectorD computeHii(TMatrixD A)
virtual void printLog(Option_t *option="")
void exportModelFortran(std::ofstream *sourcefile) const
Export the model in Fortran langage in a file.
virtual ~TLinearRegression()
Constructor with a dataserver and a PMML File.
TMatrixD getParametersCovariance() const
Return the Q2 statistic.
Definition TLinearRegression.h:116
Double_t _dr2a
The Adjusted R2 quality.
Definition TLinearRegression.h:60
void exportModelPython(std::ofstream *sourcefile) const
Export the model in Python langage in a file (not yet implemented)
Definition TLinearRegression.h:155
Definition TModeler.h:63
ROOT.
Definition TAnisp.h:164