ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MethodSVM.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Marcin Wolter, Andrzej Zemla
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : MethodSVM *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation *
12  * *
13  * Authors (alphabetical): *
14  * Marcin Wolter <Marcin.Wolter@cern.ch> - IFJ PAN, Krakow, Poland *
15  * Andrzej Zemla <azemla@cern.ch> - IFJ PAN, Krakow, Poland *
16  * (IFJ PAN: Henryk Niewodniczanski Inst. Nucl. Physics, Krakow, Poland) *
17  * *
18  * Introduction of regression by: *
19  * Krzysztof Danielowski <danielow@cern.ch> - IFJ PAN & AGH, Krakow, Poland *
20  * Kamil Kraszewski <kalq@cern.ch> - IFJ PAN & UJ, Krakow, Poland *
21  * Maciej Kruk <mkruk@cern.ch> - IFJ PAN & AGH, Krakow, Poland *
22  * *
23  * *
24  * Copyright (c) 2005: *
25  * CERN, Switzerland *
26  * MPI-K Heidelberg, Germany *
27  * PAN, Krakow, Poland *
28  * *
29  * Redistribution and use in source and binary forms, with or without *
30  * modification, are permitted according to the terms listed in LICENSE *
31  * (http://tmva.sourceforge.net/LICENSE) *
32  **********************************************************************************/
33 
34 //_______________________________________________________________________
35 //
36 // SMO Platt's SVM classifier with Keerthi & Shavade improvements
37 //_______________________________________________________________________
38 
39 #include "Riostream.h"
40 #include "TFile.h"
41 #include "TVectorD.h"
42 #include "TMath.h"
43 
44 #include "TMVA/ClassifierFactory.h"
45 #ifndef ROOT_TMVA_MethodSVM
46 #include "TMVA/MethodSVM.h"
47 #endif
48 #ifndef ROOT_TMVA_Tools
49 #include "TMVA/Tools.h"
50 #endif
51 #ifndef ROOT_TMVA_Timer
52 #include "TMVA/Timer.h"
53 #endif
54 
55 #ifndef ROOT_TMVA_SVWorkingSet
56 #include "TMVA/SVWorkingSet.h"
57 #endif
58 
59 #ifndef ROOT_TMVA_SVEvent
60 #include "TMVA/SVEvent.h"
61 #endif
62 
63 #ifndef ROOT_TMVA_SVKernelFunction
64 #include "TMVA/SVKernelFunction.h"
65 #endif
66 
67 #include "TMVA/DataSet.h"
68 #include "TMVA/DataSetInfo.h"
69 #include "TMVA/Event.h"
70 #include "TMVA/MethodBase.h"
71 #include "TMVA/MsgLogger.h"
72 #include "TMVA/Types.h"
73 
74 #include <string>
75 
76 using std::vector;
77 
78 //const Int_t basketsize__ = 1280000;
79 REGISTER_METHOD(SVM)
80 
81 ClassImp(TMVA::MethodSVM)
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// standard constructor
85 
86 TMVA::MethodSVM::MethodSVM( const TString& jobName, const TString& methodTitle, DataSetInfo& theData,
87  const TString& theOption, TDirectory* theTargetDir )
88  : MethodBase( jobName, Types::kSVM, methodTitle, theData, theOption, theTargetDir )
89  , fCost(0)
90  , fTolerance(0)
91  , fMaxIter(0)
92  , fNSubSets(0)
93  , fBparm(0)
94  , fGamma(0)
95  , fWgSet(0)
96  , fInputData(0)
97  , fSupportVectors(0)
98  , fSVKernelFunction(0)
99  , fMinVars(0)
100  , fMaxVars(0)
101  , fDoubleSigmaSquared(0)
102  , fOrder(0)
103  , fTheta(0)
104  , fKappa(0)
105 {
106 }
107 
108 ////////////////////////////////////////////////////////////////////////////////
109 /// constructor from weight file
110 
111 TMVA::MethodSVM::MethodSVM( DataSetInfo& theData, const TString& theWeightFile, TDirectory* theTargetDir )
112  : MethodBase( Types::kSVM, theData, theWeightFile, theTargetDir )
113  , fCost(0)
114  , fTolerance(0)
115  , fMaxIter(0)
116  , fNSubSets(0)
117  , fBparm(0)
118  , fGamma(0)
119  , fWgSet(0)
120  , fInputData(0)
121  , fSupportVectors(0)
122  , fSVKernelFunction(0)
123  , fMinVars(0)
124  , fMaxVars(0)
125  , fDoubleSigmaSquared(0)
126  , fOrder(0)
127  , fTheta(0)
128  , fKappa(0)
129 {
130 }
131 
132 ////////////////////////////////////////////////////////////////////////////////
133 /// destructor
134 
136 {
137  if (fInputData !=0) { delete fInputData; fInputData=0; }
138  if (fSupportVectors !=0 ) { delete fSupportVectors; fSupportVectors = 0; }
139  if (fWgSet !=0) { delete fWgSet; fWgSet=0; }
140  if (fSVKernelFunction !=0 ) { delete fSVKernelFunction; fSVKernelFunction = 0; }
141 }
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 /// SVM can handle classification with 2 classes and regression with one regression-target
145 
147 {
148  if (type == Types::kClassification && numberClasses == 2) return kTRUE;
149  if (type == Types::kRegression && numberTargets == 1) return kTRUE;
150  return kFALSE;
151 }
152 
153 ////////////////////////////////////////////////////////////////////////////////
154 /// default initialisation
155 
157 {
158  // SVM always uses normalised input variables
159  SetNormalised( kTRUE );
160 
161  // Helge: do not book a event vector of given size but rather fill the vector
162  // later with pus_back. Anyway, this is NOT what is time consuming in
163  // SVM and it allows to skip totally events with weights == 0 ;)
164  fInputData = new std::vector<TMVA::SVEvent*>(0);
165  fSupportVectors = new std::vector<TMVA::SVEvent*>(0);
166 }
167 
168 ////////////////////////////////////////////////////////////////////////////////
169 /// declare options available for this method
170 
172 {
173  // for gaussian kernel parameter(s)
174  DeclareOptionRef( fGamma = 1., "Gamma", "RBF kernel parameter: Gamma (size of the Kernel)");
175 
176  DeclareOptionRef( fCost, "C", "Cost parameter" );
177  if (DoRegression()) {
178  fCost = 0.002;
179  }else{
180  fCost = 1.;
181  }
182  DeclareOptionRef( fTolerance = 0.01, "Tol", "Tolerance parameter" ); //should be fixed
183  DeclareOptionRef( fMaxIter = 1000, "MaxIter", "Maximum number of training loops" );
184 
185 }
186 
187 ////////////////////////////////////////////////////////////////////////////////
188 /// options that are used ONLY for the READER to ensure backward compatibility
189 
191 {
193  DeclareOptionRef( fNSubSets = 1, "NSubSets", "Number of training subsets" );
194  DeclareOptionRef( fTheKernel = "Gauss", "Kernel", "Uses kernel function");
195  // for gaussian kernel parameter(s)
196  DeclareOptionRef( fDoubleSigmaSquared = 2., "Sigma", "Kernel parameter: sigma");
197  // for polynomiarl kernel parameter(s)
198  DeclareOptionRef( fOrder = 3, "Order", "Polynomial Kernel parameter: polynomial order");
199  // for sigmoid kernel parameters
200  DeclareOptionRef( fTheta = 1., "Theta", "Sigmoid Kernel parameter: theta");
201  DeclareOptionRef( fKappa = 1., "Kappa", "Sigmoid Kernel parameter: kappa");
202 }
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 /// option post processing (if necessary)
206 
208 {
209  if (IgnoreEventsWithNegWeightsInTraining()) {
210  Log() << kFATAL << "Mechanism to ignore events with negative weights in training not yet available for method: "
211  << GetMethodTypeName()
212  << " --> please remove \"IgnoreNegWeightsInTraining\" option from booking string."
213  << Endl;
214  }
215 }
216 
217 ////////////////////////////////////////////////////////////////////////////////
218 /// Train SVM
219 
221 {
222  Data()->SetCurrentType(Types::kTraining);
223 
224  Log() << kDEBUG << "Create event vector"<< Endl;
225  for (Int_t ievt=0; ievt<Data()->GetNEvents(); ievt++){
226  if (GetEvent(ievt)->GetWeight() != 0)
227  fInputData->push_back(new SVEvent(GetEvent(ievt), fCost, DataInfo().IsSignal(GetEvent(ievt))));
228  }
229 
230  fSVKernelFunction = new SVKernelFunction(fGamma);
231 
232  Log()<< kINFO << "Building SVM Working Set...with "<<fInputData->size()<<" event instances"<< Endl;
233  Timer bldwstime( GetName());
234  fWgSet = new SVWorkingSet( fInputData, fSVKernelFunction,fTolerance, DoRegression() );
235  Log() << kINFO <<"Elapsed time for Working Set build: "<< bldwstime.GetElapsedTime()<<Endl;
236 
237  // timing
238  Timer timer( GetName() );
239  Log() << kINFO << "Sorry, no computing time forecast available for SVM, please wait ..." << Endl;
240 
241  fWgSet->Train(fMaxIter);
242 
243  Log() << kINFO << "Elapsed time: " << timer.GetElapsedTime()
244  << " " << Endl;
245 
246  fBparm = fWgSet->GetBpar();
247  fSupportVectors = fWgSet->GetSupportVectors();
248 
249 
250  delete fWgSet;
251  fWgSet=0;
252 
253  // for (UInt_t i=0; i<fInputData->size();i++) delete fInputData->at(i);
254  delete fInputData;
255  fInputData=0;
256 }
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 /// write configuration to xml file
260 
261 void TMVA::MethodSVM::AddWeightsXMLTo( void* parent ) const
262 {
263  void* wght = gTools().AddChild(parent, "Weights");
264  gTools().AddAttr(wght,"fBparm",fBparm);
265  gTools().AddAttr(wght,"fGamma",fGamma);
266  gTools().AddAttr(wght,"NSupVec",fSupportVectors->size());
267 
268  for (std::vector<TMVA::SVEvent*>::iterator veciter=fSupportVectors->begin();
269  veciter!=fSupportVectors->end() ; ++veciter ) {
270  TVectorD temp(GetNvar()+4);
271  temp[0] = (*veciter)->GetNs();
272  temp[1] = (*veciter)->GetTypeFlag();
273  temp[2] = (*veciter)->GetAlpha();
274  temp[3] = (*veciter)->GetAlpha_p();
275  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++)
276  temp[ivar+4] = (*(*veciter)->GetDataVector())[ivar];
277  gTools().WriteTVectorDToXML(wght,"SupportVector",&temp);
278  }
279  // write max/min data values
280  void* maxnode = gTools().AddChild(wght, "Maxima");
281  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++)
282  gTools().AddAttr(maxnode, "Var"+gTools().StringFromInt(ivar), GetXmax(ivar));
283  void* minnode = gTools().AddChild(wght, "Minima");
284  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++)
285  gTools().AddAttr(minnode, "Var"+gTools().StringFromInt(ivar), GetXmin(ivar));
286 }
287 
288 ////////////////////////////////////////////////////////////////////////////////
289 
291 {
292  gTools().ReadAttr( wghtnode, "fBparm",fBparm );
293  gTools().ReadAttr( wghtnode, "fGamma",fGamma);
294  UInt_t fNsupv=0;
295  gTools().ReadAttr( wghtnode, "NSupVec",fNsupv );
296 
297  Float_t alpha=0.;
298  Float_t alpha_p = 0.;
299 
300  Int_t typeFlag=-1;
301  // UInt_t ns = 0;
302  std::vector<Float_t>* svector = new std::vector<Float_t>(GetNvar());
303 
304  if (fMaxVars!=0) delete fMaxVars;
305  fMaxVars = new TVectorD( GetNvar() );
306  if (fMinVars!=0) delete fMinVars;
307  fMinVars = new TVectorD( GetNvar() );
308  if (fSupportVectors!=0) {
309  for (vector< SVEvent* >::iterator it = fSupportVectors->begin(); it!=fSupportVectors->end(); ++it)
310  delete *it;
311  delete fSupportVectors;
312  }
313  fSupportVectors = new std::vector<TMVA::SVEvent*>(0);
314  void* supportvectornode = gTools().GetChild(wghtnode);
315  for (UInt_t ievt = 0; ievt < fNsupv; ievt++) {
316  TVectorD temp(GetNvar()+4);
317  gTools().ReadTVectorDFromXML(supportvectornode,"SupportVector",&temp);
318  // ns=(UInt_t)temp[0];
319  typeFlag=(int)temp[1];
320  alpha=temp[2];
321  alpha_p=temp[3];
322  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++) (*svector)[ivar]=temp[ivar+4];
323 
324  fSupportVectors->push_back(new SVEvent(svector,alpha,alpha_p,typeFlag));
325  supportvectornode = gTools().GetNextChild(supportvectornode);
326  }
327 
328  void* maxminnode = supportvectornode;
329  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++)
330  gTools().ReadAttr( maxminnode,"Var"+gTools().StringFromInt(ivar),(*fMaxVars)[ivar]);
331  maxminnode = gTools().GetNextChild(maxminnode);
332  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++)
333  gTools().ReadAttr( maxminnode,"Var"+gTools().StringFromInt(ivar),(*fMinVars)[ivar]);
334  if (fSVKernelFunction!=0) delete fSVKernelFunction;
335  fSVKernelFunction = new SVKernelFunction(fGamma);
336  delete svector;
337 }
338 
339 ////////////////////////////////////////////////////////////////////////////////
340 ///TODO write IT
341 /// write training sample (TTree) to file
342 
344 {
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 
349 void TMVA::MethodSVM::ReadWeightsFromStream( std::istream& istr )
350 {
351  if (fSupportVectors !=0) { delete fSupportVectors; fSupportVectors = 0;}
352  fSupportVectors = new std::vector<TMVA::SVEvent*>(0);
353 
354  // read configuration from input stream
355  istr >> fBparm;
356 
357  UInt_t fNsupv;
358  // coverity[tainted_data_argument]
359  istr >> fNsupv;
360  fSupportVectors->reserve(fNsupv);
361 
362  Float_t typeTalpha=0.;
363  Float_t alpha=0.;
364  Int_t typeFlag=-1;
365  UInt_t ns = 0;
366  std::vector<Float_t>* svector = new std::vector<Float_t>(GetNvar());
367 
368  fMaxVars = new TVectorD( GetNvar() );
369  fMinVars = new TVectorD( GetNvar() );
370 
371  for (UInt_t ievt = 0; ievt < fNsupv; ievt++) {
372  istr>>ns;
373  istr>>typeTalpha;
374  typeFlag = typeTalpha<0?-1:1;
375  alpha = typeTalpha<0?-typeTalpha:typeTalpha;
376  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++) istr >> svector->at(ivar);
377 
378  fSupportVectors->push_back(new SVEvent(svector,alpha,typeFlag,ns));
379  }
380 
381  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++) istr >> (*fMaxVars)[ivar];
382 
383  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++) istr >> (*fMinVars)[ivar];
384 
385  delete fSVKernelFunction;
386  if (fTheKernel == "Gauss" ) {
387  fSVKernelFunction = new SVKernelFunction(1/fDoubleSigmaSquared);
388  }
389  else {
391  if(fTheKernel == "Linear") k = SVKernelFunction::kLinear;
392  else if (fTheKernel == "Polynomial") k = SVKernelFunction::kPolynomial;
393  else if (fTheKernel == "Sigmoid" ) k = SVKernelFunction::kSigmoidal;
394  else {
395  Log() << kFATAL <<"Unknown kernel function found in weight file!" << Endl;
396  }
397  fSVKernelFunction = new SVKernelFunction();
398  fSVKernelFunction->setCompatibilityParams(k, fOrder, fTheta, fKappa);
399  }
400  delete svector;
401 }
402 
403 ////////////////////////////////////////////////////////////////////////////////
404 /// TODO write IT
405 
407 {
408 }
409 
410 ////////////////////////////////////////////////////////////////////////////////
411 /// returns MVA value for given event
412 
414 {
415  Double_t myMVA = 0;
416 
417  // TODO: avoid creation of a new SVEvent every time (Joerg)
418  SVEvent* ev = new SVEvent( GetEvent(), 0. ); // check for specificators
419 
420  for (UInt_t ievt = 0; ievt < fSupportVectors->size() ; ievt++) {
421  myMVA += ( fSupportVectors->at(ievt)->GetAlpha()
422  * fSupportVectors->at(ievt)->GetTypeFlag()
423  * fSVKernelFunction->Evaluate( fSupportVectors->at(ievt), ev ) );
424  }
425 
426  delete ev;
427 
428  myMVA -= fBparm;
429 
430  // cannot determine error
431  NoErrorCalc(err, errUpper);
432 
433  // 08/12/09: changed sign here to make results agree with convention signal=1
434  return 1.0/(1.0 + TMath::Exp(myMVA));
435 }
436 ////////////////////////////////////////////////////////////////////////////////
437 
438 const std::vector<Float_t>& TMVA::MethodSVM::GetRegressionValues()
439 {
440  if( fRegressionReturnVal == NULL )
441  fRegressionReturnVal = new std::vector<Float_t>();
442  fRegressionReturnVal->clear();
443 
444  Double_t myMVA = 0;
445 
446  const Event *baseev = GetEvent();
447  SVEvent* ev = new SVEvent( baseev,0. ); //check for specificators
448 
449  for (UInt_t ievt = 0; ievt < fSupportVectors->size() ; ievt++) {
450  myMVA += ( fSupportVectors->at(ievt)->GetDeltaAlpha()
451  *fSVKernelFunction->Evaluate( fSupportVectors->at(ievt), ev ) );
452  }
453  myMVA += fBparm;
454  Event * evT = new Event(*baseev);
455  evT->SetTarget(0,myMVA);
456 
457  const Event* evT2 = GetTransformationHandler().InverseTransform( evT );
458 
459  fRegressionReturnVal->push_back(evT2->GetTarget(0));
460 
461  delete evT;
462 
463  delete ev;
464 
465  return *fRegressionReturnVal;
466 }
467 
468 ////////////////////////////////////////////////////////////////////////////////
469 /// write specific classifier response
470 
471 void TMVA::MethodSVM::MakeClassSpecific( std::ostream& fout, const TString& className ) const
472 {
473  const int fNsupv = fSupportVectors->size();
474  fout << " // not implemented for class: \"" << className << "\"" << std::endl;
475  fout << " float fBparameter;" << std::endl;
476  fout << " int fNOfSuppVec;" << std::endl;
477  fout << " static float fAllSuppVectors[][" << fNsupv << "];" << std::endl;
478  fout << " static float fAlphaTypeCoef[" << fNsupv << "];" << std::endl;
479  fout << std::endl;
480  fout << " // Kernel parameter(s) " << std::endl;
481  fout << " float fGamma;" << std::endl;
482  fout << "};" << std::endl;
483  fout << "" << std::endl;
484 
485  //Initialize function definition
486  fout << "inline void " << className << "::Initialize() " << std::endl;
487  fout << "{" << std::endl;
488  fout << " fBparameter = " << fBparm << ";" << std::endl;
489  fout << " fNOfSuppVec = " << fNsupv << ";" << std::endl;
490  fout << " fGamma = " << fGamma << ";" <<std::endl;
491  fout << "}" << std::endl;
492  fout << std::endl;
493 
494  // GetMvaValue__ function defninition
495  fout << "inline double " << className << "::GetMvaValue__(const std::vector<double>& inputValues ) const" << std::endl;
496  fout << "{" << std::endl;
497  fout << " double mvaval = 0; " << std::endl;
498  fout << " double temp = 0; " << std::endl;
499  fout << std::endl;
500  fout << " for (int ievt = 0; ievt < fNOfSuppVec; ievt++ ){" << std::endl;
501  fout << " temp = 0;" << std::endl;
502  fout << " for ( unsigned int ivar = 0; ivar < GetNvar(); ivar++ ) {" << std::endl;
503 
504  fout << " temp += (fAllSuppVectors[ivar][ievt] - inputValues[ivar]) " << std::endl;
505  fout << " * (fAllSuppVectors[ivar][ievt] - inputValues[ivar]); " << std::endl;
506  fout << " }" << std::endl;
507  fout << " mvaval += fAlphaTypeCoef[ievt] * exp( -fGamma * temp ); " << std::endl;
508 
509  fout << " }" << std::endl;
510  fout << " mvaval -= fBparameter;" << std::endl;
511  fout << " return 1./(1. + exp(mvaval));" << std::endl;
512  fout << "}" << std::endl;
513  fout << "// Clean up" << std::endl;
514  fout << "inline void " << className << "::Clear() " << std::endl;
515  fout << "{" << std::endl;
516  fout << " // nothing to clear " << std::endl;
517  fout << "}" << std::endl;
518  fout << "" << std::endl;
519 
520  // define support vectors
521  fout << "float " << className << "::fAlphaTypeCoef[] =" << std::endl;
522  fout << "{ ";
523  for (Int_t isv = 0; isv < fNsupv; isv++) {
524  fout << fSupportVectors->at(isv)->GetDeltaAlpha() * fSupportVectors->at(isv)->GetTypeFlag();
525  if (isv < fNsupv-1) fout << ", ";
526  }
527  fout << " };" << std::endl << std::endl;
528 
529  fout << "float " << className << "::fAllSuppVectors[][" << fNsupv << "] =" << std::endl;
530  fout << "{";
531  for (UInt_t ivar = 0; ivar < GetNvar(); ivar++) {
532  fout << std::endl;
533  fout << " { ";
534  for (Int_t isv = 0; isv < fNsupv; isv++){
535  fout << fSupportVectors->at(isv)->GetDataVector()->at(ivar);
536  if (isv < fNsupv-1) fout << ", ";
537  }
538  fout << " }";
539  if (ivar < GetNvar()-1) fout << ", " << std::endl;
540  else fout << std::endl;
541  }
542  fout << "};" << std::endl<< std::endl;
543 }
544 
545 ////////////////////////////////////////////////////////////////////////////////
546 /// get help message text
547 ///
548 /// typical length of text line:
549 /// "|--------------------------------------------------------------|"
550 
552 {
553  Log() << Endl;
554  Log() << gTools().Color("bold") << "--- Short description:" << gTools().Color("reset") << Endl;
555  Log() << Endl;
556  Log() << "The Support Vector Machine (SVM) builds a hyperplance separating" << Endl;
557  Log() << "signal and background events (vectors) using the minimal subset of " << Endl;
558  Log() << "all vectors used for training (support vectors). The extension to" << Endl;
559  Log() << "the non-linear case is performed by mapping input vectors into a " << Endl;
560  Log() << "higher-dimensional feature space in which linear separation is " << Endl;
561  Log() << "possible. The use of the kernel functions thereby eliminates the " << Endl;
562  Log() << "explicit transformation to the feature space. The implemented SVM " << Endl;
563  Log() << "algorithm performs the classification tasks using linear, polynomial, " << Endl;
564  Log() << "Gaussian and sigmoidal kernel functions. The Gaussian kernel allows " << Endl;
565  Log() << "to apply any discriminant shape in the input space." << Endl;
566  Log() << Endl;
567  Log() << gTools().Color("bold") << "--- Performance optimisation:" << gTools().Color("reset") << Endl;
568  Log() << Endl;
569  Log() << "SVM is a general purpose non-linear classification method, which " << Endl;
570  Log() << "does not require data preprocessing like decorrelation or Principal " << Endl;
571  Log() << "Component Analysis. It generalises quite well and can handle analyses " << Endl;
572  Log() << "with large numbers of input variables." << Endl;
573  Log() << Endl;
574  Log() << gTools().Color("bold") << "--- Performance tuning via configuration options:" << gTools().Color("reset") << Endl;
575  Log() << Endl;
576  Log() << "Optimal performance requires primarily a proper choice of the kernel " << Endl;
577  Log() << "parameters (the width \"Sigma\" in case of Gaussian kernel) and the" << Endl;
578  Log() << "cost parameter \"C\". The user must optimise them empirically by running" << Endl;
579  Log() << "SVM several times with different parameter sets. The time needed for " << Endl;
580  Log() << "each evaluation scales like the square of the number of training " << Endl;
581  Log() << "events so that a coarse preliminary tuning should be performed on " << Endl;
582  Log() << "reduced data sets." << Endl;
583 }
void Train(void)
Train SVM.
Definition: MethodSVM.cxx:220
ClassImp(TMVA::MethodSVM) TMVA
standard constructor
Definition: MethodSVM.cxx:81
void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility
Definition: MethodSVM.cxx:190
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
float Float_t
Definition: RtypesCore.h:53
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
EAnalysisType
Definition: Types.h:124
Basic string class.
Definition: TString.h:137
void MakeClassSpecific(std::ostream &, const TString &) const
write specific classifier response
Definition: MethodSVM.cxx:471
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
void SetTarget(UInt_t itgt, Float_t value)
set the target value (dimension itgt) to value
Definition: Event.cxx:354
void ProcessOptions()
option post processing (if necessary)
Definition: MethodSVM.cxx:207
void AddAttr(void *node, const char *, const T &value, Int_t precision=16)
Definition: Tools.h:308
void * AddChild(void *parent, const char *childname, const char *content=0, bool isRootNode=false)
add child node
Definition: Tools.cxx:1134
Tools & gTools()
Definition: Tools.cxx:79
TStopwatch timer
Definition: pirndm.C:37
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
std::vector< std::vector< double > > Data
void DeclareOptions()
declare options available for this method
Definition: MethodSVM.cxx:171
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:24
TString GetElapsedTime(Bool_t Scientific=kTRUE)
Definition: Timer.cxx:131
void AddWeightsXMLTo(void *parent) const
write configuration to xml file
Definition: MethodSVM.cxx:261
void ReadWeightsFromXML(void *wghtnode)
Definition: MethodSVM.cxx:290
void ReadWeightsFromStream(std::istream &istr)
Definition: MethodSVM.cxx:349
void Init(void)
default initialisation
Definition: MethodSVM.cxx:156
unsigned int UInt_t
Definition: RtypesCore.h:42
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets)
SVM can handle classification with 2 classes and regression with one regression-target.
Definition: MethodSVM.cxx:146
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
returns MVA value for given event
Definition: MethodSVM.cxx:413
MethodSVM(const TString &jobName, const TString &methodTitle, DataSetInfo &theData, const TString &theOption="", TDirectory *theTargetDir=0)
Double_t Exp(Double_t x)
Definition: TMath.h:495
double Double_t
Definition: RtypesCore.h:55
Describe directory structure in memory.
Definition: TDirectory.h:44
int type
Definition: TGX11.cxx:120
void WriteTVectorDToXML(void *node, const char *name, TVectorD *vec)
Definition: Tools.cxx:1267
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
const TString & Color(const TString &)
human readable color strings
Definition: Tools.cxx:837
Float_t GetTarget(UInt_t itgt) const
Definition: Event.h:101
#define REGISTER_METHOD(CLASS)
for example
virtual ~MethodSVM(void)
destructor
Definition: MethodSVM.cxx:135
virtual void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility they are hence without any...
Definition: MethodBase.cxx:606
#define NULL
Definition: Rtypes.h:82
void WriteWeightsToStream(TFile &fout) const
TODO write IT write training sample (TTree) to file.
Definition: MethodSVM.cxx:343
void ReadTVectorDFromXML(void *node, const char *name, TVectorD *vec)
Definition: Tools.cxx:1275
const Bool_t kTRUE
Definition: Rtypes.h:91
void GetHelpMessage() const
get help message text
Definition: MethodSVM.cxx:551
Definition: math.cpp:60
const std::vector< Float_t > & GetRegressionValues()
Definition: MethodSVM.cxx:438