ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
VariablePCATransform.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Joerg Stelzer, Helge Voss, Eckhard von Toerne
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : VariablePCATransform *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation (see header for description) *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Peter Speckmayer <Peter.Speckmayer@cern.ch> - CERN, Switzerland *
16  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - CERN, Switzerland *
17  * Eckhard v. Toerne <evt@uni-bonn.de> - U of Bonn, Germany *
18  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
19  * *
20  * Copyright (c) 2005-2011: *
21  * CERN, Switzerland *
22  * MPI-K Heidelberg, Germany *
23  * U. of Bonn, Germany *
24  * *
25  * Redistribution and use in source and binary forms, with or without *
26  * modification, are permitted according to the terms listed in LICENSE *
27  * (http://tmva.sourceforge.net/LICENSE) *
28  **********************************************************************************/
29 
31 
32 #include "TMVA/DataSet.h"
33 #include "TMVA/Event.h"
34 #include "TMVA/MsgLogger.h"
35 #include "TMVA/Tools.h"
36 #include "TMVA/Types.h"
37 
38 #include "TMatrixD.h"
39 #include "TMatrixDBase.h"
40 #include "TPrincipal.h"
41 #include "TVectorD.h"
42 #include "TVectorF.h"
43 
44 #include <iostream>
45 #include <iomanip>
46 #include <stdexcept>
47 #include <algorithm>
48 
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 /// constructor
53 
54 TMVA::VariablePCATransform::VariablePCATransform( DataSetInfo& dsi )
55 : VariableTransformBase( dsi, Types::kPCA, "PCA" )
56 {
57 }
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 /// destructor
61 
63 {
64  for (UInt_t i=0; i<fMeanValues.size(); i++) {
65  if (fMeanValues.at(i) != 0) delete fMeanValues.at(i);
66  if (fEigenVectors.at(i) != 0) delete fEigenVectors.at(i);
67  }
68 }
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 /// initialization of the transformation.
72 /// Has to be called in the preparation and not in the constructor,
73 /// since the number of classes it not known at construction, but
74 /// only after the creation of the DataSet which might be later.
75 
77 {
78 }
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// calculate the principal components using the ROOT class TPrincipal
82 /// and the normalization
83 
84 Bool_t TMVA::VariablePCATransform::PrepareTransformation (const std::vector<Event*>& events)
85 {
86  Initialize();
87 
88  if (!IsEnabled() || IsCreated()) return kTRUE;
89 
90  Log() << kINFO << "Preparing the Principle Component (PCA) transformation..." << Endl;
91 
92  UInt_t inputSize = fGet.size();
93 
94  SetNVariables(inputSize);
95 
96  // TPrincipal doesn't support PCA transformation for 1 or less variables
97  if (inputSize <= 1) {
98  Log() << kFATAL << "Cannot perform PCA transformation for " << inputSize << " variable only" << Endl;
99  return kFALSE;
100  }
101 
102  if (inputSize > 200) {
103  Log() << kINFO << "----------------------------------------------------------------------------"
104  << Endl;
105  Log() << kINFO
106  << ": More than 200 variables, will not calculate PCA!" << Endl;
107  Log() << kINFO << "----------------------------------------------------------------------------"
108  << Endl;
109  return kFALSE;
110  }
111 
112  CalculatePrincipalComponents( events );
113 
114  SetCreated( kTRUE );
115 
116  return kTRUE;
117 }
118 
119 ////////////////////////////////////////////////////////////////////////////////
120 /// apply the principal component analysis
121 
122 const TMVA::Event* TMVA::VariablePCATransform::Transform( const Event* const ev, Int_t cls ) const
123 {
124  if (!IsCreated()) return 0;
125 
126 // const Int_t inputSize = fGet.size();
127 // const UInt_t nCls = GetNClasses();
128 
129  // if we have more than one class, take the last PCA analysis where all classes are combined if
130  // the cls parameter is outside the defined classes
131  // If there is only one class, then no extra class for all events of all classes has to be created
132 
133  //if (cls < 0 || cls > GetNClasses()) cls = (fMeanValues.size()==1?0:2);//( GetNClasses() == 1 ? 0 : 1 ); ;
134  // EVT this is a workaround to address the reader problem with transforma and EvaluateMVA(std::vector<float/double> ,...)
135  if (cls < 0 || cls >= (int) fMeanValues.size()) cls = fMeanValues.size()-1;
136  // EVT workaround end
137 
138  // Perform PCA and put it into PCAed events tree
139 
140  if (fTransformedEvent==0 ) {
141  fTransformedEvent = new Event();
142  }
143 
144  std::vector<Float_t> input;
145  std::vector<Char_t> mask;
146  std::vector<Float_t> principalComponents;
147 
148  Bool_t hasMaskedEntries = GetInput( ev, input, mask );
149 
150  if( hasMaskedEntries ){ // targets might be masked (for events where the targets have not been computed yet)
151  UInt_t numMasked = std::count(mask.begin(), mask.end(), (Char_t)kTRUE);
152  UInt_t numOK = std::count(mask.begin(), mask.end(), (Char_t)kFALSE);
153  if( numMasked>0 && numOK>0 ){
154  Log() << kFATAL << "You mixed variables and targets in the decorrelation transformation. This is not possible." << Endl;
155  }
156  SetOutput( fTransformedEvent, input, mask, ev );
157  return fTransformedEvent;
158  }
159 
160  X2P( principalComponents, input, cls );
161  SetOutput( fTransformedEvent, principalComponents, mask, ev );
162 
163  return fTransformedEvent;
164 }
165 
166 ////////////////////////////////////////////////////////////////////////////////
167 /// apply the principal component analysis
168 /// TODO: implementation of inverse transformation
169 /// Log() << kFATAL << "Inverse transformation for PCA transformation not yet implemented. Hence, this transformation cannot be applied together with regression. Please contact the authors if necessary." << Endl;
170 
172 {
173  if (!IsCreated()) return 0;
174 // const Int_t inputSize = fGet.size();
175  const UInt_t nCls = GetNClasses();
176  //UInt_t evCls = ev->GetClass();
177 
178  // if we have more than one class, take the last PCA analysis where all classes are combined if
179  // the cls parameter is outside the defined classes
180  // If there is only one class, then no extra class for all events of all classes has to be created
181  if (cls < 0 || UInt_t(cls) > nCls) cls = (fMeanValues.size()==1?0:2);//( GetNClasses() == 1 ? 0 : 1 ); ;
182  // Perform PCA and put it into PCAed events tree
183 
184  if (fBackTransformedEvent==0 ) fBackTransformedEvent = new Event();
185 
186  std::vector<Float_t> principalComponents;
187  std::vector<Char_t> mask;
188  std::vector<Float_t> output;
189 
190  GetInput( ev, principalComponents, mask, kTRUE );
191  P2X( output, principalComponents, cls );
192  SetOutput( fBackTransformedEvent, output, mask, ev, kTRUE );
193 
194  return fBackTransformedEvent;
195 }
196 
197 ////////////////////////////////////////////////////////////////////////////////
198 /// calculate the principal components for the signal and the background data
199 /// it uses the MakePrincipal method of ROOT's TPrincipal class
200 
201 void TMVA::VariablePCATransform::CalculatePrincipalComponents( const std::vector< Event*>& events )
202 {
203  UInt_t nvars = 0, ntgts = 0, nspcts = 0;
204  CountVariableTypes( nvars, ntgts, nspcts );
205  if( nvars>0 && ntgts>0 )
206  Log() << kFATAL << "Variables and targets cannot be mixed in PCA transformation." << Endl;
207 
208  const Int_t inputSize = fGet.size();
209 
210  // if we have more than one class, add another PCA analysis which combines all classes
211  const UInt_t nCls = GetNClasses();
212  const UInt_t maxPCA = (nCls<=1) ? nCls : nCls+1;
213 
214  // PCA [signal/background/class x/class y/... /all classes]
215  std::vector<TPrincipal*> pca(maxPCA);
216  for (UInt_t i=0; i<maxPCA; i++) pca[i] = new TPrincipal(nvars,"");
217 
218  // !! Not normalizing and not storing input data, for performance reasons. Should perhaps restore normalization.
219  // But this can be done afterwards by adding a normalisation transformation (user defined)
220 
221  Long64_t ievt, entries = events.size();
222  Double_t *dvec = new Double_t[inputSize];
223 
224  std::vector<Float_t> input;
225  std::vector<Char_t> mask;
226  for (ievt=0; ievt<entries; ievt++) {
227  const Event* ev = events[ievt];
228  UInt_t cls = ev->GetClass();
229 
230  Bool_t hasMaskedEntries = GetInput( ev, input, mask );
231  if (hasMaskedEntries){
232  Log() << kWARNING << "Print event which triggers an error" << Endl;
233  ev->Print(Log());
234  Log() << kFATAL << "Masked entries found in event read in when calculating the principal components for the PCA transformation." << Endl;
235  }
236 
237  UInt_t iinp = 0;
238  for( std::vector<Float_t>::iterator itInp = input.begin(), itInpEnd = input.end(); itInp != itInpEnd; ++itInp )
239  {
240  Float_t value = (*itInp);
241  dvec[iinp] = (Double_t)value;
242  ++iinp;
243  }
244 
245  pca.at(cls)->AddRow( dvec );
246  if (nCls > 1) pca.at(maxPCA-1)->AddRow( dvec );
247  }
248 
249  // delete possible leftovers
250  for (UInt_t i=0; i<fMeanValues.size(); i++) if (fMeanValues[i] != 0) delete fMeanValues[i];
251  for (UInt_t i=0; i<fEigenVectors.size(); i++) if (fEigenVectors[i] != 0) delete fEigenVectors[i];
252  fMeanValues.resize(maxPCA,0);
253  fEigenVectors.resize(maxPCA,0);
254 
255  for (UInt_t i=0; i<maxPCA; i++ ) {
256  pca.at(i)->MakePrincipals();
257 
258  // retrieve mean values, eigenvectors and sigmas
259  fMeanValues[i] = new TVectorD( *(pca.at(i)->GetMeanValues()) ); // need to copy since we want to own
260  fEigenVectors[i] = new TMatrixD( *(pca.at(i)->GetEigenVectors()) );
261  }
262 
263  for (UInt_t i=0; i<maxPCA; i++) delete pca.at(i);
264  delete [] dvec;
265 }
266 
267 ////////////////////////////////////////////////////////////////////////////////
268 /// Calculate the principal components from the original data vector
269 /// x, and return it in p (function extracted from TPrincipal::X2P)
270 /// It's the users responsibility to make sure that both x and p are
271 /// of the right size (i.e., memory must be allocated for p)
272 
273 void TMVA::VariablePCATransform::X2P( std::vector<Float_t>& pc, const std::vector<Float_t>& x, Int_t cls ) const
274 {
275  const Int_t nInput = x.size();
276  pc.assign(nInput,0);
277 
278  for (Int_t i = 0; i < nInput; i++) {
279  Double_t pv = 0;
280  for (Int_t j = 0; j < nInput; j++)
281  pv += (((Double_t)x.at(j)) - (*fMeanValues.at(cls))(j)) * (*fEigenVectors.at(cls))(j,i);
282  pc[i] = pv;
283  }
284 }
285 
286 ////////////////////////////////////////////////////////////////////////////////
287 /// Perform the back-transformation from the principal components
288 /// pc, and return x
289 /// It's the users responsibility to make sure that both x and pc are
290 /// of the right size (i.e., memory must be allocated for p)
291 
292 void TMVA::VariablePCATransform::P2X( std::vector<Float_t>& x, const std::vector<Float_t>& pc, Int_t cls ) const
293 {
294  const Int_t nInput = pc.size();
295  x.assign(nInput,0);
296 
297  for (Int_t i = 0; i < nInput; i++) {
298  Double_t xv = 0;
299  for (Int_t j = 0; j < nInput; j++)
300  xv += (((Double_t)pc.at(j)) * (*fEigenVectors.at(cls))(i,j) ) + (*fMeanValues.at(cls))(j);
301  x[i] = xv;
302  }
303 }
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// write mean values to stream
307 
309 {
310  for (Int_t sbType=0; sbType<2; sbType++) {
311  o << "# PCA mean values " << std::endl;
312  const TVectorD* means = fMeanValues[sbType];
313  o << (sbType==0 ? "Signal" : "Background") << " " << means->GetNrows() << std::endl;
314  for (Int_t row = 0; row<means->GetNrows(); row++) {
315  o << std::setprecision(12) << std::setw(20) << (*means)[row];
316  }
317  o << std::endl;
318  }
319  o << "##" << std::endl;
320 
321  // write eigenvectors to stream
322  for (Int_t sbType=0; sbType<2; sbType++) {
323  o << "# PCA eigenvectors " << std::endl;
324  const TMatrixD* mat = fEigenVectors[sbType];
325  o << (sbType==0 ? "Signal" : "Background") << " " << mat->GetNrows() << " x " << mat->GetNcols() << std::endl;
326  for (Int_t row = 0; row<mat->GetNrows(); row++) {
327  for (Int_t col = 0; col<mat->GetNcols(); col++) {
328  o << std::setprecision(12) << std::setw(20) << (*mat)[row][col] << " ";
329  }
330  o << std::endl;
331  }
332  }
333  o << "##" << std::endl;
334 }
335 
336 ////////////////////////////////////////////////////////////////////////////////
337 /// create XML description of PCA transformation
338 
340  void* trfxml = gTools().AddChild(parent, "Transform");
341  gTools().AddAttr(trfxml, "Name", "PCA");
342 
344 
345  // write mean values to stream
346  for (UInt_t sbType=0; sbType<fMeanValues.size(); sbType++) {
347  void* meanxml = gTools().AddChild( trfxml, "Statistics");
348  const TVectorD* means = fMeanValues[sbType];
349  gTools().AddAttr( meanxml, "Class", (sbType==0 ? "Signal" :(sbType==1 ? "Background":"Combined")) );
350  gTools().AddAttr( meanxml, "ClassIndex", sbType );
351  gTools().AddAttr( meanxml, "NRows", means->GetNrows() );
352  TString meansdef = "";
353  for (Int_t row = 0; row<means->GetNrows(); row++)
354  meansdef += gTools().StringFromDouble((*means)[row]) + " ";
355  gTools().AddRawLine( meanxml, meansdef );
356  }
357 
358  // write eigenvectors to stream
359  for (UInt_t sbType=0; sbType<fEigenVectors.size(); sbType++) {
360  void* evxml = gTools().AddChild( trfxml, "Eigenvectors");
361  const TMatrixD* mat = fEigenVectors[sbType];
362  gTools().AddAttr( evxml, "Class", (sbType==0 ? "Signal" :(sbType==1 ? "Background":"Combined") ) );
363  gTools().AddAttr( evxml, "ClassIndex", sbType );
364  gTools().AddAttr( evxml, "NRows", mat->GetNrows() );
365  gTools().AddAttr( evxml, "NCols", mat->GetNcols() );
366  TString evdef = "";
367  for (Int_t row = 0; row<mat->GetNrows(); row++)
368  for (Int_t col = 0; col<mat->GetNcols(); col++)
369  evdef += gTools().StringFromDouble((*mat)[row][col]) + " ";
370  gTools().AddRawLine( evxml, evdef );
371  }
372 }
373 
374 ////////////////////////////////////////////////////////////////////////////////
375 /// Read the transformation matrices from the xml node
376 
378 {
379  Int_t nrows, ncols;
380  UInt_t clsIdx;
381  TString classtype;
382  TString nodeName;
383 
384  Bool_t newFormat = kFALSE;
385 
386  void* inpnode = NULL;
387 
388  inpnode = gTools().GetChild(trfnode, "Selection"); // new xml format
389  if( inpnode!=NULL )
390  newFormat = kTRUE; // new xml format
391 
392  if( newFormat ){
393  // ------------- new format --------------------
394  // read input
396 
397  }
398 
399  void* ch = gTools().GetChild(trfnode);
400  while (ch) {
401  nodeName = gTools().GetName(ch);
402  if (nodeName == "Statistics") {
403  // read mean values
404  gTools().ReadAttr(ch, "Class", classtype);
405  gTools().ReadAttr(ch, "ClassIndex", clsIdx);
406  gTools().ReadAttr(ch, "NRows", nrows);
407 
408  // set the correct size
409  if (fMeanValues.size()<=clsIdx) fMeanValues.resize(clsIdx+1,0);
410  if (fMeanValues[clsIdx]==0) fMeanValues[clsIdx] = new TVectorD( nrows );
411  fMeanValues[clsIdx]->ResizeTo( nrows );
412 
413  // now read vector entries
414  std::stringstream s(gTools().GetContent(ch));
415  for (Int_t row = 0; row<nrows; row++) s >> (*fMeanValues[clsIdx])(row);
416  }
417  else if ( nodeName == "Eigenvectors" ) {
418  // Read eigenvectors
419  gTools().ReadAttr(ch, "Class", classtype);
420  gTools().ReadAttr(ch, "ClassIndex", clsIdx);
421  gTools().ReadAttr(ch, "NRows", nrows);
422  gTools().ReadAttr(ch, "NCols", ncols);
423 
424  if (fEigenVectors.size()<=clsIdx) fEigenVectors.resize(clsIdx+1,0);
425  if (fEigenVectors[clsIdx]==0) fEigenVectors[clsIdx] = new TMatrixD( nrows, ncols );
426  fEigenVectors[clsIdx]->ResizeTo( nrows, ncols );
427 
428  // now read matrix entries
429  std::stringstream s(gTools().GetContent(ch));
430  for (Int_t row = 0; row<nrows; row++)
431  for (Int_t col = 0; col<ncols; col++)
432  s >> (*fEigenVectors[clsIdx])[row][col];
433  } // done reading eigenvectors
434  ch = gTools().GetNextChild(ch);
435  }
436 
437  SetCreated();
438 }
439 
440 ////////////////////////////////////////////////////////////////////////////////
441 /// Read mean values from input stream
442 
443 void TMVA::VariablePCATransform::ReadTransformationFromStream( std::istream& istr, const TString& classname )
444 {
445  char buf[512];
446  istr.getline(buf,512);
447  TString strvar, dummy;
448  Int_t nrows(0), ncols(0);
449  UInt_t classIdx=(classname=="signal"?0:1);
450 
451  for (UInt_t i=0; i<fMeanValues.size(); i++) {
452  if (fMeanValues.at(i) != 0) delete fMeanValues.at(i);
453  if (fEigenVectors.at(i) != 0) delete fEigenVectors.at(i);
454  }
455  fMeanValues.resize(3);
456  fEigenVectors.resize(3);
457 
458  Log() << kINFO << "VariablePCATransform::ReadTransformationFromStream(): " << Endl;
459 
460  while (!(buf[0]=='#'&& buf[1]=='#')) { // if line starts with ## return
461  char* p = buf;
462  while (*p==' ' || *p=='\t') p++; // 'remove' leading whitespace
463  if (*p=='#' || *p=='\0') {
464  istr.getline(buf,512);
465  continue; // if comment or empty line, read the next line
466  }
467  std::stringstream sstr(buf);
468  sstr >> strvar;
469  if (strvar=="signal" || strvar=="background") {
470 
471  sstr >> nrows;
472  Int_t sbType = (strvar=="signal" ? 0 : 1);
473 
474  if (fMeanValues[sbType] == 0) fMeanValues[sbType] = new TVectorD( nrows );
475  else fMeanValues[sbType]->ResizeTo( nrows );
476 
477  // now read vector entries
478  for (Int_t row = 0; row<nrows; row++) istr >> (*fMeanValues[sbType])(row);
479 
480  } // done reading vector
481 
482  istr.getline(buf,512); // reading the next line
483  }
484 
485  // Read eigenvectors from input stream
486  istr.getline(buf,512);
487  while (!(buf[0]=='#'&& buf[1]=='#')) { // if line starts with ## return
488  char* p = buf;
489  while(*p==' ' || *p=='\t') p++; // 'remove' leading whitespace
490  if (*p=='#' || *p=='\0') {
491  istr.getline(buf,512);
492  continue; // if comment or empty line, read the next line
493  }
494  std::stringstream sstr(buf);
495  sstr >> strvar;
496  if (strvar=="signal" || strvar=="background") {
497 
498  // coverity[tainted_data_argument]
499  sstr >> nrows >> dummy >> ncols;
500  Int_t sbType = (strvar=="signal" ? 0 : 1);
501 
502  if (fEigenVectors[sbType] == 0) fEigenVectors[sbType] = new TMatrixD( nrows, ncols );
503  else fEigenVectors[sbType]->ResizeTo( nrows, ncols );
504 
505  // now read matrix entries
506  for (Int_t row = 0; row<fEigenVectors[sbType]->GetNrows(); row++) {
507  for (Int_t col = 0; col<fEigenVectors[sbType]->GetNcols(); col++) {
508  istr >> (*fEigenVectors[sbType])[row][col];
509  }
510  }
511 
512  } // done reading matrix
513  istr.getline(buf,512); // reading the next line
514  }
515  fMeanValues[2] = new TVectorD( *fMeanValues[classIdx] );
516  fEigenVectors[2] = new TMatrixD( *fEigenVectors[classIdx] );
517 
518  SetCreated();
519 }
520 
521 ////////////////////////////////////////////////////////////////////////////////
522 /// creates C++ code fragment of the PCA transform for inclusion in standalone C++ class
523 
524 void TMVA::VariablePCATransform::MakeFunction( std::ostream& fout, const TString& fcncName,
525  Int_t part, UInt_t trCounter, Int_t )
526 {
527  UInt_t nvar = fEigenVectors[0]->GetNrows();
528 
529  // creates a PCA transformation function
530  UInt_t numC = fMeanValues.size();
531  if (part==1) {
532  fout << std::endl;
533  fout << " void X2P_"<<trCounter<<"( const double*, double*, int ) const;" << std::endl;
534  fout << " double fMeanValues_"<<trCounter<<"["<<numC<<"]["
535  << fMeanValues[0]->GetNrows() << "];" << std::endl; // mean values
536  fout << " double fEigenVectors_"<<trCounter<<"["<<numC<<"]["
537  << fEigenVectors[0]->GetNrows() << "]["
538  << fEigenVectors[0]->GetNcols() <<"];" << std::endl; // eigenvectors
539  fout << std::endl;
540  }
541 
542  // sanity check
543  if (numC>1){
544  if (fMeanValues[0]->GetNrows() != fMeanValues[1]->GetNrows() ||
545  fEigenVectors[0]->GetNrows() != fEigenVectors[1]->GetNrows() ||
546  fEigenVectors[0]->GetNcols() != fEigenVectors[1]->GetNcols()) {
547  Log() << kFATAL << "<MakeFunction> Mismatch in vector/matrix dimensions" << Endl;
548  }
549  }
550 
551  if (part==2) {
552 
553  fout << std::endl;
554  fout << "//_______________________________________________________________________" << std::endl;
555  fout << "inline void " << fcncName << "::X2P_"<<trCounter<<"( const double* x, double* p, int index ) const" << std::endl;
556  fout << "{" << std::endl;
557  fout << " // Calculate the principal components from the original data vector" << std::endl;
558  fout << " // x, and return it in p (function extracted from TPrincipal::X2P)" << std::endl;
559  fout << " // It's the users responsibility to make sure that both x and p are" << std::endl;
560  fout << " // of the right size (i.e., memory must be allocated for p)." << std::endl;
561  fout << " const int nVar = " << nvar << ";" << std::endl;
562  fout << std::endl;
563  fout << " for (int i = 0; i < nVar; i++) {" << std::endl;
564  fout << " p[i] = 0;" << std::endl;
565  fout << " for (int j = 0; j < nVar; j++) p[i] += (x[j] - fMeanValues_"<<trCounter<<"[index][j]) * fEigenVectors_"<<trCounter<<"[index][j][i];" << std::endl;
566  fout << " }" << std::endl;
567  fout << "}" << std::endl;
568  fout << std::endl;
569  fout << "//_______________________________________________________________________" << std::endl;
570  fout << "inline void " << fcncName << "::InitTransform_"<<trCounter<<"()" << std::endl;
571  fout << "{" << std::endl;
572  fout << " // PCA transformation, initialisation" << std::endl;
573 
574  // fill vector of mean values
575  fout << " // initialise vector of mean values" << std::endl;
576  std::streamsize dp = fout.precision();
577  for (UInt_t index=0; index<numC; index++) {
578  for (int i=0; i<fMeanValues[index]->GetNrows(); i++) {
579  fout << " fMeanValues_"<<trCounter<<"["<<index<<"]["<<i<<"] = " << std::setprecision(12)
580  << (*fMeanValues[index])(i) << ";" << std::endl;
581  }
582  }
583 
584  // fill matrix of eigenvectors
585  fout << std::endl;
586  fout << " // initialise matrix of eigenvectors" << std::endl;
587  for (UInt_t index=0; index<numC; index++) {
588  for (int i=0; i<fEigenVectors[index]->GetNrows(); i++) {
589  for (int j=0; j<fEigenVectors[index]->GetNcols(); j++) {
590  fout << " fEigenVectors_"<<trCounter<<"["<<index<<"]["<<i<<"]["<<j<<"] = " << std::setprecision(12)
591  << (*fEigenVectors[index])(i,j) << ";" << std::endl;
592  }
593  }
594  }
595  fout << std::setprecision(dp);
596  fout << "}" << std::endl;
597  fout << std::endl;
598  fout << "//_______________________________________________________________________" << std::endl;
599  fout << "inline void " << fcncName << "::Transform_"<<trCounter<<"( std::vector<double>& iv, int cls ) const" << std::endl;
600  fout << "{" << std::endl;
601  fout << " // PCA transformation" << std::endl;
602  fout << " const int nVar = " << nvar << ";" << std::endl;
603  fout << " double *dv = new double[nVar];" << std::endl;
604  fout << " double *rv = new double[nVar];" << std::endl;
605  fout << " if (cls < 0 || cls > "<<GetNClasses()<<") {"<< std::endl;
606  fout << " if ("<<GetNClasses()<<" > 1 ) cls = "<<GetNClasses()<<";"<< std::endl;
607  fout << " else cls = "<<(numC==1?0:2)<<";"<< std::endl;
608  fout << " }"<< std::endl;
609 
610  VariableTransformBase::MakeFunction(fout, fcncName, 0, trCounter, 0 );
611 
612  fout << " for (int ivar=0; ivar<nVar; ivar++) dv[ivar] = iv[indicesGet.at(ivar)];" << std::endl;
613 
614  fout << std::endl;
615  fout << " // Perform PCA and put it into PCAed events tree" << std::endl;
616  fout << " this->X2P_"<<trCounter<<"( dv, rv, cls );" << std::endl;
617  fout << " for (int ivar=0; ivar<nVar; ivar++) iv[indicesPut.at(ivar)] = rv[ivar];" << std::endl;
618 
619  fout << std::endl;
620  fout << " delete [] dv;" << std::endl;
621  fout << " delete [] rv;" << std::endl;
622  fout << "}" << std::endl;
623  }
624 }
tuple row
Definition: mrt.py:26
Principal Components Analysis (PCA)
Definition: TPrincipal.h:28
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
long long Long64_t
Definition: RtypesCore.h:69
Bool_t PrepareTransformation(const std::vector< Event * > &)
calculate the principal components using the ROOT class TPrincipal and the normalization ...
float Float_t
Definition: RtypesCore.h:53
void ReadTransformationFromStream(std::istream &, const TString &)
Read mean values from input stream.
virtual void MakeFunction(std::ostream &fout, const TString &fncName, Int_t part, UInt_t trCounter, Int_t cls)=0
getinput and setoutput equivalent
Int_t GetNrows() const
Definition: TVectorT.h:81
virtual void AttachXMLTo(void *parent)=0
create XML description the transformation (write out info of selected variables)
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
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
virtual void AttachXMLTo(void *parent)
create XML description of PCA transformation
void Print(std::ostream &o) const
print method
Definition: Event.cxx:346
Tools & gTools()
Definition: Tools.cxx:79
Double_t x[n]
Definition: legend1.C:17
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:24
Bool_t AddRawLine(void *node, const char *raw)
XML helpers.
Definition: Tools.cxx:1198
virtual void ReadFromXML(void *trfnode)=0
Read the input variables from the XML node.
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:24
TString StringFromDouble(Double_t d)
string tools
Definition: Tools.cxx:1241
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
virtual const Event * InverseTransform(const Event *const, Int_t cls) const
apply the principal component analysis TODO: implementation of inverse transformation Log() << kFATAL...
ClassImp(TMVA::VariablePCATransform) TMVA
constructor
unsigned int UInt_t
Definition: RtypesCore.h:42
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
virtual void ReadFromXML(void *trfnode)
Read the transformation matrices from the xml node.
double Double_t
Definition: RtypesCore.h:55
Int_t GetNcols() const
Definition: TMatrixTBase.h:137
static RooMathCoreReg dummy
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
UInt_t GetClass() const
Definition: Event.h:86
void Initialize()
initialization of the transformation.
char Char_t
Definition: RtypesCore.h:29
std::vector< TMatrixD * > fEigenVectors
virtual const Event * Transform(const Event *const, Int_t cls) const
apply the principal component analysis
#define NULL
Definition: Rtypes.h:82
const char * GetName(void *node)
XML helpers.
Definition: Tools.cxx:1190
virtual ~VariablePCATransform(void)
destructor
void P2X(std::vector< Float_t > &, const std::vector< Float_t > &, Int_t cls) const
Perform the back-transformation from the principal components pc, and return x It's the users respons...
std::vector< TVectorD * > fMeanValues
static void output(int code)
Definition: gifencode.c:226
const Bool_t kTRUE
Definition: Rtypes.h:91
void WriteTransformationToStream(std::ostream &) const
write mean values to stream
float value
Definition: math.cpp:443
Definition: math.cpp:60
void X2P(std::vector< Float_t > &, const std::vector< Float_t > &, Int_t cls) const
Calculate the principal components from the original data vector x, and return it in p (function extr...
virtual void MakeFunction(std::ostream &fout, const TString &fncName, Int_t part, UInt_t trCounter, Int_t cls)
creates C++ code fragment of the PCA transform for inclusion in standalone C++ class ...
void CalculatePrincipalComponents(const std::vector< Event * > &)
calculate the principal components for the signal and the background data it uses the MakePrincipal m...