Logo ROOT   6.08/07
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 
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 
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 
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 
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  std::ostringstream oss;
234  ev->Print(oss);
235  Log() << oss.str();
236  Log() << kFATAL << "Masked entries found in event read in when calculating the principal components for the PCA transformation." << Endl;
237  }
238 
239  UInt_t iinp = 0;
240  for( std::vector<Float_t>::iterator itInp = input.begin(), itInpEnd = input.end(); itInp != itInpEnd; ++itInp )
241  {
242  Float_t value = (*itInp);
243  dvec[iinp] = (Double_t)value;
244  ++iinp;
245  }
246 
247  pca.at(cls)->AddRow( dvec );
248  if (nCls > 1) pca.at(maxPCA-1)->AddRow( dvec );
249  }
250 
251  // delete possible leftovers
252  for (UInt_t i=0; i<fMeanValues.size(); i++) if (fMeanValues[i] != 0) delete fMeanValues[i];
253  for (UInt_t i=0; i<fEigenVectors.size(); i++) if (fEigenVectors[i] != 0) delete fEigenVectors[i];
254  fMeanValues.resize(maxPCA,0);
255  fEigenVectors.resize(maxPCA,0);
256 
257  for (UInt_t i=0; i<maxPCA; i++ ) {
258  pca.at(i)->MakePrincipals();
259 
260  // retrieve mean values, eigenvectors and sigmas
261  fMeanValues[i] = new TVectorD( *(pca.at(i)->GetMeanValues()) ); // need to copy since we want to own
262  fEigenVectors[i] = new TMatrixD( *(pca.at(i)->GetEigenVectors()) );
263  }
264 
265  for (UInt_t i=0; i<maxPCA; i++) delete pca.at(i);
266  delete [] dvec;
267 }
268 
269 ////////////////////////////////////////////////////////////////////////////////
270 /// Calculate the principal components from the original data vector
271 /// x, and return it in p (function extracted from TPrincipal::X2P)
272 /// It's the users responsibility to make sure that both x and p are
273 /// of the right size (i.e., memory must be allocated for p)
274 
275 void TMVA::VariablePCATransform::X2P( std::vector<Float_t>& pc, const std::vector<Float_t>& x, Int_t cls ) const
276 {
277  const Int_t nInput = x.size();
278  pc.assign(nInput,0);
279 
280  for (Int_t i = 0; i < nInput; i++) {
281  Double_t pv = 0;
282  for (Int_t j = 0; j < nInput; j++)
283  pv += (((Double_t)x.at(j)) - (*fMeanValues.at(cls))(j)) * (*fEigenVectors.at(cls))(j,i);
284  pc[i] = pv;
285  }
286 }
287 
288 ////////////////////////////////////////////////////////////////////////////////
289 /// Perform the back-transformation from the principal components
290 /// pc, and return x
291 /// It's the users responsibility to make sure that both x and pc are
292 /// of the right size (i.e., memory must be allocated for p)
293 
294 void TMVA::VariablePCATransform::P2X( std::vector<Float_t>& x, const std::vector<Float_t>& pc, Int_t cls ) const
295 {
296  const Int_t nInput = pc.size();
297  x.assign(nInput,0);
298 
299  for (Int_t i = 0; i < nInput; i++) {
300  Double_t xv = 0;
301  for (Int_t j = 0; j < nInput; j++)
302  xv += (((Double_t)pc.at(j)) * (*fEigenVectors.at(cls))(i,j) ) + (*fMeanValues.at(cls))(j);
303  x[i] = xv;
304  }
305 }
306 
307 ////////////////////////////////////////////////////////////////////////////////
308 /// write mean values to stream
309 
311 {
312  for (Int_t sbType=0; sbType<2; sbType++) {
313  o << "# PCA mean values " << std::endl;
314  const TVectorD* means = fMeanValues[sbType];
315  o << (sbType==0 ? "Signal" : "Background") << " " << means->GetNrows() << std::endl;
316  for (Int_t row = 0; row<means->GetNrows(); row++) {
317  o << std::setprecision(12) << std::setw(20) << (*means)[row];
318  }
319  o << std::endl;
320  }
321  o << "##" << std::endl;
322 
323  // write eigenvectors to stream
324  for (Int_t sbType=0; sbType<2; sbType++) {
325  o << "# PCA eigenvectors " << std::endl;
326  const TMatrixD* mat = fEigenVectors[sbType];
327  o << (sbType==0 ? "Signal" : "Background") << " " << mat->GetNrows() << " x " << mat->GetNcols() << std::endl;
328  for (Int_t row = 0; row<mat->GetNrows(); row++) {
329  for (Int_t col = 0; col<mat->GetNcols(); col++) {
330  o << std::setprecision(12) << std::setw(20) << (*mat)[row][col] << " ";
331  }
332  o << std::endl;
333  }
334  }
335  o << "##" << std::endl;
336 }
337 
338 ////////////////////////////////////////////////////////////////////////////////
339 /// create XML description of PCA transformation
340 
342  void* trfxml = gTools().AddChild(parent, "Transform");
343  gTools().AddAttr(trfxml, "Name", "PCA");
344 
346 
347  // write mean values to stream
348  for (UInt_t sbType=0; sbType<fMeanValues.size(); sbType++) {
349  void* meanxml = gTools().AddChild( trfxml, "Statistics");
350  const TVectorD* means = fMeanValues[sbType];
351  gTools().AddAttr( meanxml, "Class", (sbType==0 ? "Signal" :(sbType==1 ? "Background":"Combined")) );
352  gTools().AddAttr( meanxml, "ClassIndex", sbType );
353  gTools().AddAttr( meanxml, "NRows", means->GetNrows() );
354  TString meansdef = "";
355  for (Int_t row = 0; row<means->GetNrows(); row++)
356  meansdef += gTools().StringFromDouble((*means)[row]) + " ";
357  gTools().AddRawLine( meanxml, meansdef );
358  }
359 
360  // write eigenvectors to stream
361  for (UInt_t sbType=0; sbType<fEigenVectors.size(); sbType++) {
362  void* evxml = gTools().AddChild( trfxml, "Eigenvectors");
363  const TMatrixD* mat = fEigenVectors[sbType];
364  gTools().AddAttr( evxml, "Class", (sbType==0 ? "Signal" :(sbType==1 ? "Background":"Combined") ) );
365  gTools().AddAttr( evxml, "ClassIndex", sbType );
366  gTools().AddAttr( evxml, "NRows", mat->GetNrows() );
367  gTools().AddAttr( evxml, "NCols", mat->GetNcols() );
368  TString evdef = "";
369  for (Int_t row = 0; row<mat->GetNrows(); row++)
370  for (Int_t col = 0; col<mat->GetNcols(); col++)
371  evdef += gTools().StringFromDouble((*mat)[row][col]) + " ";
372  gTools().AddRawLine( evxml, evdef );
373  }
374 }
375 
376 ////////////////////////////////////////////////////////////////////////////////
377 /// Read the transformation matrices from the xml node
378 
380 {
381  Int_t nrows, ncols;
382  UInt_t clsIdx;
383  TString classtype;
384  TString nodeName;
385 
386  Bool_t newFormat = kFALSE;
387 
388  void* inpnode = NULL;
389 
390  inpnode = gTools().GetChild(trfnode, "Selection"); // new xml format
391  if( inpnode!=NULL )
392  newFormat = kTRUE; // new xml format
393 
394  if( newFormat ){
395  // ------------- new format --------------------
396  // read input
398 
399  }
400 
401  void* ch = gTools().GetChild(trfnode);
402  while (ch) {
403  nodeName = gTools().GetName(ch);
404  if (nodeName == "Statistics") {
405  // read mean values
406  gTools().ReadAttr(ch, "Class", classtype);
407  gTools().ReadAttr(ch, "ClassIndex", clsIdx);
408  gTools().ReadAttr(ch, "NRows", nrows);
409 
410  // set the correct size
411  if (fMeanValues.size()<=clsIdx) fMeanValues.resize(clsIdx+1,0);
412  if (fMeanValues[clsIdx]==0) fMeanValues[clsIdx] = new TVectorD( nrows );
413  fMeanValues[clsIdx]->ResizeTo( nrows );
414 
415  // now read vector entries
416  std::stringstream s(gTools().GetContent(ch));
417  for (Int_t row = 0; row<nrows; row++) s >> (*fMeanValues[clsIdx])(row);
418  }
419  else if ( nodeName == "Eigenvectors" ) {
420  // Read eigenvectors
421  gTools().ReadAttr(ch, "Class", classtype);
422  gTools().ReadAttr(ch, "ClassIndex", clsIdx);
423  gTools().ReadAttr(ch, "NRows", nrows);
424  gTools().ReadAttr(ch, "NCols", ncols);
425 
426  if (fEigenVectors.size()<=clsIdx) fEigenVectors.resize(clsIdx+1,0);
427  if (fEigenVectors[clsIdx]==0) fEigenVectors[clsIdx] = new TMatrixD( nrows, ncols );
428  fEigenVectors[clsIdx]->ResizeTo( nrows, ncols );
429 
430  // now read matrix entries
431  std::stringstream s(gTools().GetContent(ch));
432  for (Int_t row = 0; row<nrows; row++)
433  for (Int_t col = 0; col<ncols; col++)
434  s >> (*fEigenVectors[clsIdx])[row][col];
435  } // done reading eigenvectors
436  ch = gTools().GetNextChild(ch);
437  }
438 
439  SetCreated();
440 }
441 
442 ////////////////////////////////////////////////////////////////////////////////
443 /// Read mean values from input stream
444 
445 void TMVA::VariablePCATransform::ReadTransformationFromStream( std::istream& istr, const TString& classname )
446 {
447  char buf[512];
448  istr.getline(buf,512);
449  TString strvar, dummy;
450  Int_t nrows(0), ncols(0);
451  UInt_t classIdx=(classname=="signal"?0:1);
452 
453  for (UInt_t i=0; i<fMeanValues.size(); i++) {
454  if (fMeanValues.at(i) != 0) delete fMeanValues.at(i);
455  if (fEigenVectors.at(i) != 0) delete fEigenVectors.at(i);
456  }
457  fMeanValues.resize(3);
458  fEigenVectors.resize(3);
459 
460  Log() << kINFO << "VariablePCATransform::ReadTransformationFromStream(): " << Endl;
461 
462  while (!(buf[0]=='#'&& buf[1]=='#')) { // if line starts with ## return
463  char* p = buf;
464  while (*p==' ' || *p=='\t') p++; // 'remove' leading whitespace
465  if (*p=='#' || *p=='\0') {
466  istr.getline(buf,512);
467  continue; // if comment or empty line, read the next line
468  }
469  std::stringstream sstr(buf);
470  sstr >> strvar;
471  if (strvar=="signal" || strvar=="background") {
472 
473  sstr >> nrows;
474  Int_t sbType = (strvar=="signal" ? 0 : 1);
475 
476  if (fMeanValues[sbType] == 0) fMeanValues[sbType] = new TVectorD( nrows );
477  else fMeanValues[sbType]->ResizeTo( nrows );
478 
479  // now read vector entries
480  for (Int_t row = 0; row<nrows; row++) istr >> (*fMeanValues[sbType])(row);
481 
482  } // done reading vector
483 
484  istr.getline(buf,512); // reading the next line
485  }
486 
487  // Read eigenvectors from input stream
488  istr.getline(buf,512);
489  while (!(buf[0]=='#'&& buf[1]=='#')) { // if line starts with ## return
490  char* p = buf;
491  while(*p==' ' || *p=='\t') p++; // 'remove' leading whitespace
492  if (*p=='#' || *p=='\0') {
493  istr.getline(buf,512);
494  continue; // if comment or empty line, read the next line
495  }
496  std::stringstream sstr(buf);
497  sstr >> strvar;
498  if (strvar=="signal" || strvar=="background") {
499 
500  // coverity[tainted_data_argument]
501  sstr >> nrows >> dummy >> ncols;
502  Int_t sbType = (strvar=="signal" ? 0 : 1);
503 
504  if (fEigenVectors[sbType] == 0) fEigenVectors[sbType] = new TMatrixD( nrows, ncols );
505  else fEigenVectors[sbType]->ResizeTo( nrows, ncols );
506 
507  // now read matrix entries
508  for (Int_t row = 0; row<fEigenVectors[sbType]->GetNrows(); row++) {
509  for (Int_t col = 0; col<fEigenVectors[sbType]->GetNcols(); col++) {
510  istr >> (*fEigenVectors[sbType])[row][col];
511  }
512  }
513 
514  } // done reading matrix
515  istr.getline(buf,512); // reading the next line
516  }
517  fMeanValues[2] = new TVectorD( *fMeanValues[classIdx] );
518  fEigenVectors[2] = new TMatrixD( *fEigenVectors[classIdx] );
519 
520  SetCreated();
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 /// creates C++ code fragment of the PCA transform for inclusion in standalone C++ class
525 
526 void TMVA::VariablePCATransform::MakeFunction( std::ostream& fout, const TString& fcncName,
527  Int_t part, UInt_t trCounter, Int_t )
528 {
529  UInt_t nvar = fEigenVectors[0]->GetNrows();
530 
531  // creates a PCA transformation function
532  UInt_t numC = fMeanValues.size();
533  if (part==1) {
534  fout << std::endl;
535  fout << " void X2P_"<<trCounter<<"( const double*, double*, int ) const;" << std::endl;
536  fout << " double fMeanValues_"<<trCounter<<"["<<numC<<"]["
537  << fMeanValues[0]->GetNrows() << "];" << std::endl; // mean values
538  fout << " double fEigenVectors_"<<trCounter<<"["<<numC<<"]["
539  << fEigenVectors[0]->GetNrows() << "]["
540  << fEigenVectors[0]->GetNcols() <<"];" << std::endl; // eigenvectors
541  fout << std::endl;
542  }
543 
544  // sanity check
545  if (numC>1){
546  if (fMeanValues[0]->GetNrows() != fMeanValues[1]->GetNrows() ||
547  fEigenVectors[0]->GetNrows() != fEigenVectors[1]->GetNrows() ||
548  fEigenVectors[0]->GetNcols() != fEigenVectors[1]->GetNcols()) {
549  Log() << kFATAL << "<MakeFunction> Mismatch in vector/matrix dimensions" << Endl;
550  }
551  }
552 
553  if (part==2) {
554 
555  fout << std::endl;
556  fout << "//_______________________________________________________________________" << std::endl;
557  fout << "inline void " << fcncName << "::X2P_"<<trCounter<<"( const double* x, double* p, int index ) const" << std::endl;
558  fout << "{" << std::endl;
559  fout << " // Calculate the principal components from the original data vector" << std::endl;
560  fout << " // x, and return it in p (function extracted from TPrincipal::X2P)" << std::endl;
561  fout << " // It's the users responsibility to make sure that both x and p are" << std::endl;
562  fout << " // of the right size (i.e., memory must be allocated for p)." << std::endl;
563  fout << " const int nVar = " << nvar << ";" << std::endl;
564  fout << std::endl;
565  fout << " for (int i = 0; i < nVar; i++) {" << std::endl;
566  fout << " p[i] = 0;" << std::endl;
567  fout << " for (int j = 0; j < nVar; j++) p[i] += (x[j] - fMeanValues_"<<trCounter<<"[index][j]) * fEigenVectors_"<<trCounter<<"[index][j][i];" << std::endl;
568  fout << " }" << std::endl;
569  fout << "}" << std::endl;
570  fout << std::endl;
571  fout << "//_______________________________________________________________________" << std::endl;
572  fout << "inline void " << fcncName << "::InitTransform_"<<trCounter<<"()" << std::endl;
573  fout << "{" << std::endl;
574  fout << " // PCA transformation, initialisation" << std::endl;
575 
576  // fill vector of mean values
577  fout << " // initialise vector of mean values" << std::endl;
578  std::streamsize dp = fout.precision();
579  for (UInt_t index=0; index<numC; index++) {
580  for (int i=0; i<fMeanValues[index]->GetNrows(); i++) {
581  fout << " fMeanValues_"<<trCounter<<"["<<index<<"]["<<i<<"] = " << std::setprecision(12)
582  << (*fMeanValues[index])(i) << ";" << std::endl;
583  }
584  }
585 
586  // fill matrix of eigenvectors
587  fout << std::endl;
588  fout << " // initialise matrix of eigenvectors" << std::endl;
589  for (UInt_t index=0; index<numC; index++) {
590  for (int i=0; i<fEigenVectors[index]->GetNrows(); i++) {
591  for (int j=0; j<fEigenVectors[index]->GetNcols(); j++) {
592  fout << " fEigenVectors_"<<trCounter<<"["<<index<<"]["<<i<<"]["<<j<<"] = " << std::setprecision(12)
593  << (*fEigenVectors[index])(i,j) << ";" << std::endl;
594  }
595  }
596  }
597  fout << std::setprecision(dp);
598  fout << "}" << std::endl;
599  fout << std::endl;
600  fout << "//_______________________________________________________________________" << std::endl;
601  fout << "inline void " << fcncName << "::Transform_"<<trCounter<<"( std::vector<double>& iv, int cls ) const" << std::endl;
602  fout << "{" << std::endl;
603  fout << " // PCA transformation" << std::endl;
604  fout << " const int nVar = " << nvar << ";" << std::endl;
605  fout << " double *dv = new double[nVar];" << std::endl;
606  fout << " double *rv = new double[nVar];" << std::endl;
607  fout << " if (cls < 0 || cls > "<<GetNClasses()<<") {"<< std::endl;
608  fout << " if ("<<GetNClasses()<<" > 1 ) cls = "<<GetNClasses()<<";"<< std::endl;
609  fout << " else cls = "<<(numC==1?0:2)<<";"<< std::endl;
610  fout << " }"<< std::endl;
611 
612  VariableTransformBase::MakeFunction(fout, fcncName, 0, trCounter, 0 );
613 
614  fout << " for (int ivar=0; ivar<nVar; ivar++) dv[ivar] = iv[indicesGet.at(ivar)];" << std::endl;
615 
616  fout << std::endl;
617  fout << " // Perform PCA and put it into PCAed events tree" << std::endl;
618  fout << " this->X2P_"<<trCounter<<"( dv, rv, cls );" << std::endl;
619  fout << " for (int ivar=0; ivar<nVar; ivar++) iv[indicesPut.at(ivar)] = rv[ivar];" << std::endl;
620 
621  fout << std::endl;
622  fout << " delete [] dv;" << std::endl;
623  fout << " delete [] rv;" << std::endl;
624  fout << "}" << std::endl;
625  }
626 }
Principal Components Analysis (PCA)
Definition: TPrincipal.h:28
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
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&#39;s the users respons...
long long Long64_t
Definition: RtypesCore.h:69
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 GetNcols() const
Definition: TMatrixTBase.h:137
virtual void CountVariableTypes(UInt_t &nvars, UInt_t &ntgts, UInt_t &nspcts) const
count variables, targets and spectators
TVectorT.
Definition: TMatrixTBase.h:89
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_t GetNrows() const
Definition: TVectorT.h:81
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
void WriteTransformationToStream(std::ostream &) const
write mean values to stream
void AddAttr(void *node, const char *, const T &value, Int_t precision=16)
Definition: Tools.h:309
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
virtual Bool_t GetInput(const Event *event, std::vector< Float_t > &input, std::vector< Char_t > &mask, Bool_t backTransform=kFALSE) const
select the values from the event
Double_t x[n]
Definition: legend1.C:17
virtual void SetOutput(Event *event, std::vector< Float_t > &output, std::vector< Char_t > &mask, const Event *oldEvent=0, Bool_t backTransform=kFALSE) const
select the values from the event
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
UInt_t GetClass() const
Definition: Event.h:89
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.
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:24
void CalculatePrincipalComponents(const std::vector< Event *> &)
calculate the principal components for the signal and the background data it uses the MakePrincipal m...
TString StringFromDouble(Double_t d)
string tools
Definition: Tools.cxx:1241
Bool_t PrepareTransformation(const std::vector< Event *> &)
calculate the principal components using the ROOT class TPrincipal and the normalization ...
unsigned int UInt_t
Definition: RtypesCore.h:42
void Print(std::ostream &o) const
print method
Definition: Event.cxx:348
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:296
virtual void ReadFromXML(void *trfnode)
Read the transformation matrices from the xml node.
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
static RooMathCoreReg dummy
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
virtual const Event * Transform(const Event *const, Int_t cls) const
apply the principal component analysis
void Initialize()
initialization of the transformation.
char Char_t
Definition: RtypesCore.h:29
Abstract ClassifierFactory template that handles arbitrary types.
std::vector< TMatrixD * > fEigenVectors
#define NULL
Definition: Rtypes.h:82
const char * GetName(void *node)
XML helpers.
Definition: Tools.cxx:1190
virtual ~VariablePCATransform(void)
destructor
std::vector< TVectorD * > fMeanValues
const Bool_t kTRUE
Definition: Rtypes.h:91
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 const Event * InverseTransform(const Event *const, Int_t cls) const
apply the principal component analysis TODO: implementation of inverse transformation Log() << kFATAL...
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 ...