ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
VariableDecorrTransform.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 : VariableDecorrTransform *
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  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
18  * *
19  * Copyright (c) 2005-2011: *
20  * CERN, Switzerland *
21  * MPI-K Heidelberg, Germany *
22  * U. of Bonn, Germany *
23  * *
24  * Redistribution and use in source and binary forms, with or without *
25  * modification, are permitted according to the terms listed in LICENSE *
26  * (http://tmva.sourceforge.net/LICENSE) *
27  **********************************************************************************/
28 
30 
31 #include "TMVA/DataSet.h"
32 #include "TMVA/Event.h"
33 #include "TMVA/MsgLogger.h"
34 #include "TMVA/Tools.h"
35 #include "TMVA/Types.h"
36 #include "TMVA/VariableInfo.h"
37 
38 #include "TVectorF.h"
39 #include "TVectorD.h"
40 #include "TMatrixD.h"
41 #include "TMatrixDBase.h"
42 
43 #include <iostream>
44 #include <iomanip>
45 #include <algorithm>
46 
48 
49 ////////////////////////////////////////////////////////////////////////////////
50 /// constructor
51 
52 TMVA::VariableDecorrTransform::VariableDecorrTransform( DataSetInfo& dsi )
53  : VariableTransformBase( dsi, Types::kDecorrelated, "Deco" )
54 {
55 }
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// destructor
59 
61 {
62  for (std::vector<TMatrixD*>::iterator it = fDecorrMatrices.begin(); it != fDecorrMatrices.end(); it++) {
63  if ((*it) != 0) delete (*it);
64  }
65 }
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 /// initialisation
69 
71 {
72 }
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 /// calculate the decorrelation matrix and the normalization
76 
78 {
79  Initialize();
80 
81  if (!IsEnabled() || IsCreated()) return kTRUE;
82 
83  Log() << kINFO << "Preparing the Decorrelation transformation..." << Endl;
84 
85  Int_t inputSize = fGet.size();
86  SetNVariables(inputSize);
87 
88  if (inputSize > 200) {
89  Log() << kINFO << "----------------------------------------------------------------------------"
90  << Endl;
91  Log() << kINFO
92  << ": More than 200 variables, will not calculate decorrelation matrix "
93  << "!" << Endl;
94  Log() << kINFO << "----------------------------------------------------------------------------"
95  << Endl;
96  return kFALSE;
97  }
98 
99  CalcSQRMats( events, GetNClasses() );
100 
101  SetCreated( kTRUE );
102 
103  return kTRUE;
104 }
105 
106 ////////////////////////////////////////////////////////////////////////////////
107 /// creates string with variable transformations applied
108 
110 {
111  Int_t whichMatrix = cls;
112  // if cls (the class chosen by the user) not existing, assume that user wants to
113  // have the matrix for all classes together.
114 
115  if (cls < 0 || cls > GetNClasses()) whichMatrix = GetNClasses();
116 
117  TMatrixD* m = fDecorrMatrices.at(whichMatrix);
118  if (m == 0) {
119  if (whichMatrix == GetNClasses() )
120  Log() << kFATAL << "Transformation matrix all classes is not defined"
121  << Endl;
122  else
123  Log() << kFATAL << "Transformation matrix for class " << whichMatrix << " is not defined"
124  << Endl;
125  }
126 
127  const Int_t nvar = fGet.size();
128  std::vector<TString>* strVec = new std::vector<TString>;
129 
130  // fill vector
131  for (Int_t ivar=0; ivar<nvar; ivar++) {
132  TString str( "" );
133  for (Int_t jvar=0; jvar<nvar; jvar++) {
134  str += ((*m)(ivar,jvar) > 0) ? " + " : " - ";
135 
136  Char_t type = fGet.at(jvar).first;
137  Int_t idx = fGet.at(jvar).second;
138 
139  switch( type ) {
140  case 'v':
141  str += Form( "%10.5g*[%s]", TMath::Abs((*m)(ivar,jvar)), Variables()[idx].GetLabel().Data() );
142  break;
143  case 't':
144  str += Form( "%10.5g*[%s]", TMath::Abs((*m)(ivar,jvar)), Targets()[idx].GetLabel().Data() );
145  break;
146  case 's':
147  str += Form( "%10.5g*[%s]", TMath::Abs((*m)(ivar,jvar)), Spectators()[idx].GetLabel().Data() );
148  break;
149  default:
150  Log() << kFATAL << "VariableDecorrTransform::GetTransformationStrings : unknown type '" << type << "'." << Endl;
151  }
152  }
153  strVec->push_back( str );
154  }
155 
156  return strVec;
157 }
158 
159 ////////////////////////////////////////////////////////////////////////////////
160 /// apply the decorrelation transformation
161 
163 {
164  if (!IsCreated())
165  Log() << kFATAL << "Transformation matrix not yet created"
166  << Endl;
167 
168  Int_t whichMatrix = cls;
169  // if cls (the class chosen by the user) not existing, assume that he wants to have the matrix for all classes together.
170  // EVT this is a workaround to address the reader problem with transforma and EvaluateMVA(std::vector<float/double> ,...)
171  if (cls < 0 || cls >= (int) fDecorrMatrices.size()) whichMatrix = fDecorrMatrices.size()-1;
172  //EVT workaround end
173  //if (cls < 0 || cls > GetNClasses()) {
174  // whichMatrix = GetNClasses();
175  // if (GetNClasses() == 1 ) whichMatrix = (fDecorrMatrices.size()==1?0:2);
176  //}
177 
178  TMatrixD* m = fDecorrMatrices.at(whichMatrix);
179  if (m == 0) {
180  if (whichMatrix == GetNClasses() )
181  Log() << kFATAL << "Transformation matrix all classes is not defined"
182  << Endl;
183  else
184  Log() << kFATAL << "Transformation matrix for class " << whichMatrix << " is not defined"
185  << Endl;
186  }
187 
188  if (fTransformedEvent==0 || fTransformedEvent->GetNVariables()!=ev->GetNVariables()) {
189  if (fTransformedEvent!=0) { delete fTransformedEvent; fTransformedEvent = 0; }
190  fTransformedEvent = new Event();
191  }
192 
193  // transformation to decorrelate the variables
194  const Int_t nvar = fGet.size();
195 
196  std::vector<Float_t> input;
197  std::vector<Char_t> mask; // entries with kTRUE must not be transformed
198  Bool_t hasMaskedEntries = GetInput( ev, input, mask );
199 
200  if( hasMaskedEntries ){ // targets might be masked (for events where the targets have not been computed yet)
201  UInt_t numMasked = std::count(mask.begin(), mask.end(), (Char_t)kTRUE);
202  UInt_t numOK = std::count(mask.begin(), mask.end(), (Char_t)kFALSE);
203  if( numMasked>0 && numOK>0 ){
204  Log() << kFATAL << "You mixed variables and targets in the decorrelation transformation. This is not possible." << Endl;
205  }
206  SetOutput( fTransformedEvent, input, mask, ev );
207  return fTransformedEvent;
208  }
209 
210  TVectorD vec( nvar );
211  for (Int_t ivar=0; ivar<nvar; ivar++) vec(ivar) = input.at(ivar);
212 
213  // diagonalise variable vectors
214  vec *= *m;
215 
216  input.clear();
217  for (Int_t ivar=0; ivar<nvar; ivar++) input.push_back( vec(ivar) );
218 
219  SetOutput( fTransformedEvent, input, mask, ev );
220 
221  return fTransformedEvent;
222 }
223 
224 ////////////////////////////////////////////////////////////////////////////////
225 /// apply the inverse decorrelation transformation ...
226 /// TODO : ... build the inverse transformation
227 
229 {
230  Log() << kFATAL << "Inverse transformation for decorrelation transformation not yet implemented. Hence, this transformation cannot be applied together with regression if targets should be transformed. Please contact the authors if necessary." << Endl;
231 
232 
233  return fBackTransformedEvent;
234 }
235 
236 
237 ////////////////////////////////////////////////////////////////////////////////
238 /// compute square-root matrices for signal and background
239 
240 void TMVA::VariableDecorrTransform::CalcSQRMats( const std::vector< Event*>& events, Int_t maxCls )
241 {
242  // delete old matrices if any
243  for (std::vector<TMatrixD*>::iterator it = fDecorrMatrices.begin();
244  it != fDecorrMatrices.end(); it++)
245  if (0 != (*it) ) { delete (*it); *it=0; }
246 
247 
248  // if more than one classes, then produce one matrix for all events as well (beside the matrices for each class)
249  const UInt_t matNum = (maxCls<=1)?maxCls:maxCls+1;
250  fDecorrMatrices.resize( matNum, (TMatrixD*) 0 );
251 
252  std::vector<TMatrixDSym*>* covMat = gTools().CalcCovarianceMatrices( events, maxCls, this );
253 
254 
255  for (UInt_t cls=0; cls<matNum; cls++) {
256  TMatrixD* sqrMat = gTools().GetSQRootMatrix( covMat->at(cls) );
257  if ( sqrMat==0 )
258  Log() << kFATAL << "<GetSQRMats> Zero pointer returned for SQR matrix" << Endl;
259  fDecorrMatrices[cls] = sqrMat;
260  delete (*covMat)[cls];
261  }
262  delete covMat;
263 }
264 
265 ////////////////////////////////////////////////////////////////////////////////
266 /// write the decorrelation matrix to the stream
267 
269 {
270  Int_t cls = 0;
271  Int_t dp = o.precision();
272  for (std::vector<TMatrixD*>::const_iterator itm = fDecorrMatrices.begin(); itm != fDecorrMatrices.end(); itm++) {
273  o << "# correlation matrix " << std::endl;
274  TMatrixD* mat = (*itm);
275  o << cls << " " << mat->GetNrows() << " x " << mat->GetNcols() << std::endl;
276  for (Int_t row = 0; row<mat->GetNrows(); row++) {
277  for (Int_t col = 0; col<mat->GetNcols(); col++) {
278  o << std::setprecision(12) << std::setw(20) << (*mat)[row][col] << " ";
279  }
280  o << std::endl;
281  }
282  cls++;
283  }
284  o << "##" << std::endl;
285  o << std::setprecision(dp);
286 }
287 
288 ////////////////////////////////////////////////////////////////////////////////
289 /// node attachment to parent
290 
292 {
293  void* trf = gTools().AddChild(parent, "Transform");
294  gTools().AddAttr(trf,"Name", "Decorrelation");
295 
297 
298  for (std::vector<TMatrixD*>::const_iterator itm = fDecorrMatrices.begin(); itm != fDecorrMatrices.end(); itm++) {
299  TMatrixD* mat = (*itm);
300  /*void* decmat = gTools().xmlengine().NewChild(trf, 0, "Matrix");
301  gTools().xmlengine().NewAttr(decmat,0,"Rows", gTools().StringFromInt(mat->GetNrows()) );
302  gTools().xmlengine().NewAttr(decmat,0,"Columns", gTools().StringFromInt(mat->GetNcols()) );
303 
304  std::stringstream s;
305  for (Int_t row = 0; row<mat->GetNrows(); row++) {
306  for (Int_t col = 0; col<mat->GetNcols(); col++) {
307  s << (*mat)[row][col] << " ";
308  }
309  }
310  gTools().xmlengine().AddRawLine( decmat, s.str().c_str() );*/
311  gTools().WriteTMatrixDToXML(trf,"Matrix",mat);
312  }
313 }
314 
315 ////////////////////////////////////////////////////////////////////////////////
316 /// Read the transformation matrices from the xml node
317 
319 {
320  // first delete the old matrices
321  for( std::vector<TMatrixD*>::iterator it = fDecorrMatrices.begin(); it != fDecorrMatrices.end(); it++ )
322  if( (*it) != 0 ) delete (*it);
323  fDecorrMatrices.clear();
324 
325  Bool_t newFormat = kFALSE;
326 
327  void* inpnode = NULL;
328 
329  inpnode = gTools().GetChild(trfnode, "Selection"); // new xml format
330  if( inpnode!=NULL )
331  newFormat = kTRUE; // new xml format
332 
333  void* ch = NULL;
334  if( newFormat ){
335  // ------------- new format --------------------
336  // read input
338 
339  ch = gTools().GetNextChild(inpnode);
340  }else
341  ch = gTools().GetChild(trfnode);
342 
343  // Read the transformation matrices from the xml node
344  while(ch!=0) {
345  Int_t nrows, ncols;
346  gTools().ReadAttr(ch, "Rows", nrows);
347  gTools().ReadAttr(ch, "Columns", ncols);
348  TMatrixD* mat = new TMatrixD(nrows,ncols);
349  const char* content = gTools().GetContent(ch);
350  std::stringstream s(content);
351  for (Int_t row = 0; row<nrows; row++) {
352  for (Int_t col = 0; col<ncols; col++) {
353  s >> (*mat)[row][col];
354  }
355  }
356  fDecorrMatrices.push_back(mat);
357  ch = gTools().GetNextChild(ch);
358  }
359  SetCreated();
360 }
361 
362 
363 ////////////////////////////////////////////////////////////////////////////////
364 /// Read the decorellation matrix from an input stream
365 
366 void TMVA::VariableDecorrTransform::ReadTransformationFromStream( std::istream& istr, const TString& classname )
367 {
368  char buf[512];
369  istr.getline(buf,512);
370  TString strvar, dummy;
371  Int_t nrows(0), ncols(0);
372  UInt_t classIdx=0;
373  while (!(buf[0]=='#'&& buf[1]=='#')) { // if line starts with ## return
374  char* p = buf;
375  while (*p==' ' || *p=='\t') p++; // 'remove' leading whitespace
376  if (*p=='#' || *p=='\0') {
377  istr.getline(buf,512);
378  continue; // if comment or empty line, read the next line
379  }
380  std::stringstream sstr(buf);
381 
382  sstr >> strvar;
383  if (strvar=="signal" || strvar=="background") {
384  UInt_t cls=0;
385  if(strvar=="background") cls=1;
386  if(strvar==classname) classIdx = cls;
387  // coverity[tainted_data_argument]
388  sstr >> nrows >> dummy >> ncols;
389  if (fDecorrMatrices.size() <= cls ) fDecorrMatrices.resize(cls+1);
390  if (fDecorrMatrices.at(cls) != 0) delete fDecorrMatrices.at(cls);
391  TMatrixD* mat = fDecorrMatrices.at(cls) = new TMatrixD(nrows,ncols);
392  // now read all matrix parameters
393  for (Int_t row = 0; row<mat->GetNrows(); row++) {
394  for (Int_t col = 0; col<mat->GetNcols(); col++) {
395  istr >> (*mat)[row][col];
396  }
397  }
398  } // done reading a matrix
399  istr.getline(buf,512); // reading the next line
400  }
401 
402  fDecorrMatrices.push_back( new TMatrixD(*fDecorrMatrices[classIdx]) );
403 
404  SetCreated();
405 }
406 
407 ////////////////////////////////////////////////////////////////////////////////
408 /// prints the transformation matrix
409 
411 {
412  Int_t cls = 0;
413  for (std::vector<TMatrixD*>::iterator itm = fDecorrMatrices.begin(); itm != fDecorrMatrices.end(); itm++) {
414  Log() << kINFO << "Transformation matrix "<< cls <<":" << Endl;
415  (*itm)->Print();
416  }
417 }
418 
419 ////////////////////////////////////////////////////////////////////////////////
420 /// creates C++ code fragment of the decorrelation transform for inclusion in standalone C++ class
421 
422 void TMVA::VariableDecorrTransform::MakeFunction( std::ostream& fout, const TString& fcncName, Int_t part, UInt_t trCounter, Int_t )
423 {
424  Int_t dp = fout.precision();
425 
426  UInt_t numC = fDecorrMatrices.size();
427  // creates a decorrelation function
428  if (part==1) {
429  TMatrixD* mat = fDecorrMatrices.at(0); // ToDo check if all Decorr matrices have identical dimensions
430  fout << std::endl;
431  fout << " double fDecTF_"<<trCounter<<"["<<numC<<"]["<<mat->GetNrows()<<"]["<<mat->GetNcols()<<"];" << std::endl;
432  }
433 
434  if (part==2) {
435  fout << std::endl;
436  fout << "//_______________________________________________________________________" << std::endl;
437  fout << "inline void " << fcncName << "::InitTransform_"<<trCounter<<"()" << std::endl;
438  fout << "{" << std::endl;
439  fout << " // Decorrelation transformation, initialisation" << std::endl;
440  for (UInt_t icls = 0; icls < numC; icls++){
441  TMatrixD* matx = fDecorrMatrices.at(icls);
442  for (int i=0; i<matx->GetNrows(); i++) {
443  for (int j=0; j<matx->GetNcols(); j++) {
444  fout << " fDecTF_"<<trCounter<<"["<<icls<<"]["<<i<<"]["<<j<<"] = " << std::setprecision(12) << (*matx)[i][j] << ";" << std::endl;
445  }
446  }
447  }
448  fout << "}" << std::endl;
449  fout << std::endl;
450  TMatrixD* matx = fDecorrMatrices.at(0); // ToDo check if all Decorr matrices have identicla dimensions
451  fout << "//_______________________________________________________________________" << std::endl;
452  fout << "inline void " << fcncName << "::Transform_"<<trCounter<<"( std::vector<double>& iv, int cls) const" << std::endl;
453  fout << "{" << std::endl;
454  fout << " // Decorrelation transformation" << std::endl;
455  fout << " if (cls < 0 || cls > "<<GetNClasses()<<") {"<< std::endl;
456  fout << " if ("<<GetNClasses()<<" > 1 ) cls = "<<GetNClasses()<<";"<< std::endl;
457  fout << " else cls = "<<(fDecorrMatrices.size()==1?0:2)<<";"<< std::endl;
458  fout << " }"<< std::endl;
459 
460  VariableTransformBase::MakeFunction(fout, fcncName, 0, trCounter, 0 );
461 
462  fout << " std::vector<double> tv;" << std::endl;
463  fout << " for (int i=0; i<"<<matx->GetNrows()<<";i++) {" << std::endl;
464  fout << " double v = 0;" << std::endl;
465  fout << " for (int j=0; j<"<<matx->GetNcols()<<"; j++)" << std::endl;
466  fout << " v += iv[indicesGet.at(j)] * fDecTF_"<<trCounter<<"[cls][i][j];" << std::endl;
467  fout << " tv.push_back(v);" << std::endl;
468  fout << " }" << std::endl;
469  fout << " for (int i=0; i<"<<matx->GetNrows()<<";i++) iv[indicesPut.at(i)] = tv[i];" << std::endl;
470  fout << "}" << std::endl;
471  }
472 
473  fout << std::setprecision(dp);
474 }
tuple row
Definition: mrt.py:26
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
virtual void MakeFunction(std::ostream &fout, const TString &fncName, Int_t part, UInt_t trCounter, Int_t cls)=0
getinput and setoutput equivalent
std::vector< TMatrixDSym * > * CalcCovarianceMatrices(const std::vector< Event * > &events, Int_t maxCls, VariableTransformBase *transformBase=0)
compute covariance matrices
Definition: Tools.cxx:1521
Bool_t PrepareTransformation(const std::vector< Event * > &)
calculate the decorrelation matrix and the normalization
ClassImp(TMVA::VariableDecorrTransform) TMVA
constructor
virtual void AttachXMLTo(void *parent)=0
create XML description the transformation (write out info of selected variables)
virtual void PrintTransformation(std::ostream &o)
prints the transformation matrix
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
virtual ~VariableDecorrTransform(void)
destructor
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 MakeFunction(std::ostream &fout, const TString &fncName, Int_t part, UInt_t trCounter, Int_t cls)
creates C++ code fragment of the decorrelation transform for inclusion in standalone C++ class ...
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Tools & gTools()
Definition: Tools.cxx:79
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
virtual void ReadFromXML(void *trfnode)
Read the transformation matrices from the xml node.
std::vector< std::vector< double > > Data
UInt_t GetNVariables() const
accessor to the number of variables
Definition: Event.cxx:303
virtual void ReadFromXML(void *trfnode)=0
Read the input variables from the XML node.
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:594
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:24
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
void ReadTransformationFromStream(std::istream &, const TString &)
Read the decorellation matrix from an input stream.
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
char * Form(const char *fmt,...)
const char * GetContent(void *node)
XML helpers.
Definition: Tools.cxx:1182
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
void CalcSQRMats(const std::vector< Event * > &, Int_t maxCls)
Decorrelation matrix [class0/class1/.../all classes].
void WriteTransformationToStream(std::ostream &) const
write the decorrelation matrix to the stream
std::vector< TString > * GetTransformationStrings(Int_t cls) const
creates string with variable transformations applied
virtual void AttachXMLTo(void *parent)
node attachment to parent
virtual const Event * InverseTransform(const Event *const, Int_t cls) const
apply the inverse decorrelation transformation ...
int type
Definition: TGX11.cxx:120
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
void WriteTMatrixDToXML(void *node, const char *name, TMatrixD *mat)
XML helpers.
Definition: Tools.cxx:1251
char Char_t
Definition: RtypesCore.h:29
virtual const Event * Transform(const Event *const, Int_t cls) const
apply the decorrelation transformation
#define NULL
Definition: Rtypes.h:82
TMatrixD * GetSQRootMatrix(TMatrixDSym *symMat)
square-root of symmetric matrix of course the resulting sqrtMat is also symmetric, but it's easier to treat it as a general matrix
Definition: Tools.cxx:284
const Bool_t kTRUE
Definition: Rtypes.h:91
std::vector< TMatrixD * > fDecorrMatrices
Definition: math.cpp:60