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