ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TransformationHandler.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Peter Speckmayer, Joerg Stelzer, Helge Voss, Eckhard von Toerne, Jan Therhaag
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : TransformationHandler *
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 <speckmay@mail.cern.ch> - CERN, Switzerland *
16  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - CERN, Switzerland *
17  * Jan Therhaag <Jan.Therhaag@cern.ch> - U of Bonn, Germany *
18  * Eckhard v. Toerne <evt@uni-bonn.de> - U of Bonn, Germany *
19  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
20  * *
21  * Copyright (c) 2005-2011: *
22  * CERN, Switzerland *
23  * MPI-K Heidelberg, Germany *
24  * U. of Bonn, Germany *
25  * *
26  * Redistribution and use in source and binary forms, with or without *
27  * modification, are permitted according to the terms listed in LICENSE *
28  * (http://tmva.sourceforge.net/LICENSE) *
29  **********************************************************************************/
30 
32 
33 #include "TMVA/Config.h"
34 #include "TMVA/DataSet.h"
35 #include "TMVA/DataSetInfo.h"
36 #include "TMVA/Event.h"
37 #include "TMVA/MsgLogger.h"
38 #include "TMVA/Ranking.h"
39 #include "TMVA/Tools.h"
40 #include "TMVA/Types.h"
43 #include "TMVA/VariableInfo.h"
49 
50 #include "TMath.h"
51 #include "TH1.h"
52 #include "TH2.h"
53 #include "TAxis.h"
54 #include "TProfile.h"
55 
56 #include <vector>
57 #include <iomanip>
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 /// constructor
61 
63  : fDataSetInfo(dsi),
64  fRootBaseDir(0),
65  fCallerName (callerName),
66  fLogger ( new MsgLogger(TString("TFHandler_" + callerName).Data(), kINFO) )
67 {
68  // produce one entry for each class and one entry for all classes. If there is only one class,
69  // produce only one entry
70  fNumC = (dsi.GetNClasses()<= 1) ? 1 : dsi.GetNClasses()+1;
71 
72  fVariableStats.resize( fNumC );
73  for (Int_t i=0; i<fNumC; i++ ) fVariableStats.at(i).resize(dsi.GetNVariables() + dsi.GetNTargets());
74 }
75 
76 ////////////////////////////////////////////////////////////////////////////////
77 /// destructor
78 
80 {
81  std::vector<Ranking*>::const_iterator it = fRanking.begin();
82  for (; it != fRanking.end(); it++) delete *it;
83 
84  fTransformations.SetOwner();
85  delete fLogger;
86 }
87 
88 ////////////////////////////////////////////////////////////////////////////////
89 
91 {
92  fCallerName = name;
93  fLogger->SetSource( TString("TFHandler_" + fCallerName).Data() );
94 }
95 
96 ////////////////////////////////////////////////////////////////////////////////
97 
99 {
100  TString tfname = trf->Log().GetName();
101  trf->Log().SetSource(TString(fCallerName+"_"+tfname+"_TF").Data());
102  fTransformations.Add(trf);
103  fTransformationsReferenceClasses.push_back( cls );
104  return trf;
105 }
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 
110 {
111  if (rms <= 0) {
112  Log() << kWARNING << "Variable \"" << Variable(ivar).GetExpression()
113  << "\" has zero or negative RMS^2 "
114  << "==> set to zero. Please check the variable content" << Endl;
115  rms = 0;
116  }
117 
118  VariableStat stat; stat.fMean = mean; stat.fRMS = rms; stat.fMin = min; stat.fMax = max;
119  fVariableStats.at(k).at(ivar) = stat;
120 }
121 
122 ////////////////////////////////////////////////////////////////////////////////
123 /// overrides the setting for all classes! (this is put in basically for the likelihood-method)
124 /// be careful with the usage this method
125 
127 {
128  for (UInt_t i = 0; i < fTransformationsReferenceClasses.size(); i++) {
129  fTransformationsReferenceClasses.at( i ) = cls;
130  }
131 }
132 
133 ////////////////////////////////////////////////////////////////////////////////
134 /// the transformation
135 
137 {
138  TListIter trIt(&fTransformations);
139  std::vector<Int_t>::const_iterator rClsIt = fTransformationsReferenceClasses.begin();
140  const Event* trEv = ev;
141  while (VariableTransformBase *trf = (VariableTransformBase*) trIt()) {
142  if (rClsIt == fTransformationsReferenceClasses.end()) Log() << kFATAL<< "invalid read in TransformationHandler::Transform " <<Endl;
143  trEv = trf->Transform(trEv, (*rClsIt) );
144  rClsIt++;
145  }
146  return trEv;
147 }
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 
151 const TMVA::Event* TMVA::TransformationHandler::InverseTransform( const Event* ev, Bool_t suppressIfNoTargets ) const
152 {
153  if (fTransformationsReferenceClasses.empty()){
154  //Log() << kWARNING << __FILE__ <<":InverseTransform fTransformationsReferenceClasses is empty" << Endl;
155  return ev;
156  }
157  // the inverse transformation
158  TListIter trIt(&fTransformations, kIterBackward);
159  std::vector< Int_t >::const_iterator rClsIt = fTransformationsReferenceClasses.end();
160  rClsIt--;
161  const Event* trEv = ev;
162  UInt_t nvars = 0, ntgts = 0, nspcts = 0;
163  while (VariableTransformBase *trf = (VariableTransformBase*) trIt() ) { // shouldn't be the transformation called in the inverse order for the inversetransformation?????
164  if (trf->IsCreated()) {
165  trf->CountVariableTypes( nvars, ntgts, nspcts );
166  if( !(suppressIfNoTargets && ntgts==0) )
167  trEv = trf->InverseTransform(ev, (*rClsIt) );
168  }
169  else break;
170  --rClsIt;
171  }
172  return trEv;
173 
174 
175  // TListIter trIt(&fTransformations);
176  // std::vector< Int_t >::const_iterator rClsIt = fTransformationsReferenceClasses.begin();
177  // const Event* trEv = ev;
178  // UInt_t nvars = 0, ntgts = 0, nspcts = 0;
179  // while (VariableTransformBase *trf = (VariableTransformBase*) trIt() ) { // shouldn't be the transformation called in the inverse order for the inversetransformation?????
180  // if (trf->IsCreated()) {
181  // trf->CountVariableTypes( nvars, ntgts, nspcts );
182  // if( !(suppressIfNoTargets && ntgts==0) )
183  // trEv = trf->InverseTransform(ev, (*rClsIt) );
184  // }
185  // else break;
186  // rClsIt++;
187  // }
188  // return trEv;
189 
190 }
191 
192 ////////////////////////////////////////////////////////////////////////////////
193 /// computation of transformation
194 
195 const std::vector<TMVA::Event*>* TMVA::TransformationHandler::CalcTransformations( const std::vector<Event*>& events,
196  Bool_t createNewVector )
197 {
198  if (fTransformations.GetEntries() <= 0)
199  return &events;
200 
201  // the transformedEvents are initialised with the initial events and then
202  // subsequently replaced with transformed ones. The n-th transformation will
203  // and on the events as they look like after the (n-1)-the transformation
204  // as intended for the chained transfromations
205  std::vector<Event*> *transformedEvents = new std::vector<TMVA::Event*>(events.size());
206  for ( UInt_t ievt = 0; ievt<events.size(); ievt++)
207  transformedEvents->at(ievt) = new Event(*events.at(ievt));
208 
209  TListIter trIt(&fTransformations);
210  std::vector< Int_t >::iterator rClsIt = fTransformationsReferenceClasses.begin();
211  while (VariableTransformBase *trf = (VariableTransformBase*) trIt()) {
212  if (trf->PrepareTransformation(*transformedEvents)) {
213  for (UInt_t ievt = 0; ievt<transformedEvents->size(); ievt++) { // loop through all events
214  *(*transformedEvents)[ievt] = *trf->Transform((*transformedEvents)[ievt],(*rClsIt));
215  }
216  rClsIt++;
217  }
218  }
219 
220  CalcStats(*transformedEvents);
221 
222  // plot the variables once in this transformation
223  PlotVariables(*transformedEvents);
224 
225  //sometimes, the actual transformed event vector is not used for anything but the previous
226  //CalcStat and PlotVariables calles, in that case, we delete it again (and return NULL)
227  if (!createNewVector) { // if we don't want that newly created event vector to persist, then delete it
228  for ( UInt_t ievt = 0; ievt<transformedEvents->size(); ievt++)
229  delete (*transformedEvents)[ievt];
230  delete transformedEvents;
231  transformedEvents=NULL;
232  }
233 
234  return transformedEvents; // give back the newly created event collection (containing the transformed events)
235 }
236 
237 ////////////////////////////////////////////////////////////////////////////////
238 
239 void TMVA::TransformationHandler::CalcStats (const std::vector<Event*>& events )
240 {
241  // method to calculate minimum, maximum, mean, and RMS for all
242  // variables used in the MVA
243 
244  UInt_t nevts = events.size();
245 
246  if (nevts==0)
247  Log() << kFATAL << "No events available to find min, max, mean and rms" << Endl;
248 
249  // if transformation has not been succeeded, the tree may be empty
250  const UInt_t nvar = events[0]->GetNVariables();
251  const UInt_t ntgt = events[0]->GetNTargets();
252 
253  Double_t *sumOfWeights = new Double_t[fNumC];
254  Double_t* *x2 = new Double_t*[fNumC];
255  Double_t* *x0 = new Double_t*[fNumC];
256  Double_t* *varMin = new Double_t*[fNumC];
257  Double_t* *varMax = new Double_t*[fNumC];
258 
259  for (Int_t cls=0; cls<fNumC; cls++) {
260  sumOfWeights[cls]=0;
261  x2[cls] = new Double_t[nvar+ntgt];
262  x0[cls] = new Double_t[nvar+ntgt];
263  varMin[cls] = new Double_t[nvar+ntgt];
264  varMax[cls] = new Double_t[nvar+ntgt];
265  for (UInt_t ivar=0; ivar<nvar+ntgt; ivar++) {
266  x0[cls][ivar] = x2[cls][ivar] = 0;
267  varMin[cls][ivar] = DBL_MAX;
268  varMax[cls][ivar] = -DBL_MAX;
269  }
270  }
271 
272  for (UInt_t ievt=0; ievt<nevts; ievt++) {
273  const Event* ev = events[ievt];
274  Int_t cls = ev->GetClass();
275 
276  Double_t weight = ev->GetWeight();
277  sumOfWeights[cls] += weight;
278  if (fNumC > 1 ) sumOfWeights[fNumC-1] += weight; // if more than one class, store values for all classes
279  for (UInt_t var_tgt = 0; var_tgt < 2; var_tgt++ ){ // first for variables, then for targets
280  UInt_t nloop = ( var_tgt==0?nvar:ntgt );
281  for (UInt_t ivar=0; ivar<nloop; ivar++) {
282  Double_t x = ( var_tgt==0?ev->GetValue(ivar):ev->GetTarget(ivar) );
283 
284  if (x < varMin[cls][(var_tgt*nvar)+ivar]) varMin[cls][(var_tgt*nvar)+ivar]= x;
285  if (x > varMax[cls][(var_tgt*nvar)+ivar]) varMax[cls][(var_tgt*nvar)+ivar]= x;
286 
287  x0[cls][(var_tgt*nvar)+ivar] += x*weight;
288  x2[cls][(var_tgt*nvar)+ivar] += x*x*weight;
289 
290  if (fNumC > 1) {
291  if (x < varMin[fNumC-1][(var_tgt*nvar)+ivar]) varMin[fNumC-1][(var_tgt*nvar)+ivar]= x;
292  if (x > varMax[fNumC-1][(var_tgt*nvar)+ivar]) varMax[fNumC-1][(var_tgt*nvar)+ivar]= x;
293 
294  x0[fNumC-1][(var_tgt*nvar)+ivar] += x*weight;
295  x2[fNumC-1][(var_tgt*nvar)+ivar] += x*x*weight;
296  }
297  }
298  }
299  }
300 
301 
302  // set Mean and RMS
303  for (UInt_t var_tgt = 0; var_tgt < 2; var_tgt++ ){ // first for variables, then for targets
304  UInt_t nloop = ( var_tgt==0?nvar:ntgt );
305  for (UInt_t ivar=0; ivar<nloop; ivar++) {
306  for (Int_t cls = 0; cls < fNumC; cls++) {
307  Double_t mean = x0[cls][(var_tgt*nvar)+ivar]/sumOfWeights[cls];
308  Double_t rms = TMath::Sqrt( x2[cls][(var_tgt*nvar)+ivar]/sumOfWeights[cls] - mean*mean);
309  AddStats(cls, (var_tgt*nvar)+ivar, mean, rms, varMin[cls][(var_tgt*nvar)+ivar], varMax[cls][(var_tgt*nvar)+ivar]);
310  }
311  }
312  }
313 
314  // ------ pretty output of basic statistics -------------------------------
315  // find maximum length in V (and column title)
316  UInt_t maxL = 8, maxV = 0;
317  std::vector<UInt_t> vLengths;
318  for (UInt_t ivar=0; ivar<nvar+ntgt; ivar++) {
319  if( ivar < nvar )
320  maxL = TMath::Max( (UInt_t)Variable(ivar).GetLabel().Length(), maxL );
321  else
322  maxL = TMath::Max( (UInt_t)Target(ivar-nvar).GetLabel().Length(), maxL );
323  }
324  maxV = maxL + 2;
325  // full column length
326  UInt_t clen = maxL + 4*maxV + 11;
327  for (UInt_t i=0; i<clen; i++) Log() << "-";
328  Log() << Endl;
329  // full column length
330  Log() << std::setw(maxL) << "Variable";
331  Log() << " " << std::setw(maxV) << "Mean";
332  Log() << " " << std::setw(maxV) << "RMS";
333  Log() << " " << std::setw(maxV) << "[ Min ";
334  Log() << " " << std::setw(maxV) << " Max ]" << Endl;;
335  for (UInt_t i=0; i<clen; i++) Log() << "-";
336  Log() << Endl;
337 
338  // the numbers
339  TString format = "%#11.5g";
340  for (UInt_t ivar=0; ivar<nvar+ntgt; ivar++) {
341  if( ivar < nvar )
342  Log() << std::setw(maxL) << Variable(ivar).GetLabel() << ":";
343  else
344  Log() << std::setw(maxL) << Target(ivar-nvar).GetLabel() << ":";
345  Log() << std::setw(maxV) << Form( format.Data(), GetMean(ivar) );
346  Log() << std::setw(maxV) << Form( format.Data(), GetRMS(ivar) );
347  Log() << " [" << std::setw(maxV) << Form( format.Data(), GetMin(ivar) );
348  Log() << std::setw(maxV) << Form( format.Data(), GetMax(ivar) ) << " ]";
349  Log() << Endl;
350  }
351  for (UInt_t i=0; i<clen; i++) Log() << "-";
352  Log() << Endl;
353  // ------------------------------------------------------------------------
354 
355  delete[] sumOfWeights;
356  for (Int_t cls=0; cls<fNumC; cls++) {
357  delete [] x2[cls];
358  delete [] x0[cls];
359  delete [] varMin[cls];
360  delete [] varMax[cls];
361  }
362  delete [] x2;
363  delete [] x0;
364  delete [] varMin;
365  delete [] varMax;
366 }
367 
368 ////////////////////////////////////////////////////////////////////////////////
369 /// create transformation function
370 
371 void TMVA::TransformationHandler::MakeFunction( std::ostream& fout, const TString& fncName, Int_t part ) const
372 {
373  TListIter trIt(&fTransformations);
374  std::vector< Int_t >::const_iterator rClsIt = fTransformationsReferenceClasses.begin();
375  UInt_t trCounter=1;
376  while (VariableTransformBase *trf = (VariableTransformBase*) trIt() ) {
377  trf->MakeFunction(fout, fncName, part, trCounter++, (*rClsIt) );
378  rClsIt++;
379  }
380  if (part==1) {
381  for (Int_t i=0; i<fTransformations.GetSize(); i++) {
382  fout << " void InitTransform_"<<i+1<<"();" << std::endl;
383  fout << " void Transform_"<<i+1<<"( std::vector<double> & iv, int sigOrBgd ) const;" << std::endl;
384  }
385  }
386  if (part==2) {
387  fout << std::endl;
388  fout << "//_______________________________________________________________________" << std::endl;
389  fout << "inline void " << fncName << "::InitTransform()" << std::endl;
390  fout << "{" << std::endl;
391  for (Int_t i=0; i<fTransformations.GetSize(); i++)
392  fout << " InitTransform_"<<i+1<<"();" << std::endl;
393  fout << "}" << std::endl;
394  fout << std::endl;
395  fout << "//_______________________________________________________________________" << std::endl;
396  fout << "inline void " << fncName << "::Transform( std::vector<double>& iv, int sigOrBgd ) const" << std::endl;
397  fout << "{" << std::endl;
398  for (Int_t i=0; i<fTransformations.GetSize(); i++)
399  fout << " Transform_"<<i+1<<"( iv, sigOrBgd );" << std::endl;
400 
401  fout << "}" << std::endl;
402  }
403 }
404 
405 ////////////////////////////////////////////////////////////////////////////////
406 /// return transformation name
407 
409 {
410  TString name("Id");
411  TListIter trIt(&fTransformations);
413  if ((trf = (VariableTransformBase*) trIt())) {
414  name = TString(trf->GetShortName());
415  while ((trf = (VariableTransformBase*) trIt())) name += "_" + TString(trf->GetShortName());
416  }
417  return name;
418 }
419 
420 ////////////////////////////////////////////////////////////////////////////////
421 /// incorporates transformation type into title axis (usually for histograms)
422 
424 {
425  TString xtit = info.GetTitle();
426  // indicate transformation, but not in case of single identity transform
427  if (fTransformations.GetSize() >= 1) {
428  if (fTransformations.GetSize() > 1 ||
429  ((VariableTransformBase*)GetTransformationList().Last())->GetVariableTransform() != Types::kIdentity) {
430  xtit += " (" + GetName() + ")";
431  }
432  }
433  return xtit;
434 }
435 
436 
437 ////////////////////////////////////////////////////////////////////////////////
438 /// create histograms from the input variables
439 /// - histograms for all input variables
440 /// - scatter plots for all pairs of input variables
441 
442 void TMVA::TransformationHandler::PlotVariables (const std::vector<Event*>& events, TDirectory* theDirectory )
443 {
444  if (fRootBaseDir==0 && theDirectory == 0) return;
445 
446  Log() << kINFO << "Plot event variables for ";
447  if (theDirectory !=0) Log()<< TString(theDirectory->GetName()) << Endl;
448  else Log() << GetName() << Endl;
449 
450  // extension for transformation type
451  TString transfType = "";
452  if (theDirectory == 0) {
453  transfType += "_";
454  transfType += GetName();
455  }else{ // you plot for the individual classifiers. Note, here the "statistics" still need to be calculated as you are in the testing phase
456  CalcStats(events);
457  }
458 
459  const UInt_t nvar = fDataSetInfo.GetNVariables();
460  const UInt_t ntgt = fDataSetInfo.GetNTargets();
461  const Int_t ncls = fDataSetInfo.GetNClasses();
462 
463  // Create all histograms
464  // do both, scatter and profile plots
465  std::vector<std::vector<TH1*> > hVars( ncls ); // histograms for variables
466  std::vector<std::vector<std::vector<TH2F*> > > mycorr( ncls ); // histograms for correlations
467  std::vector<std::vector<std::vector<TProfile*> > > myprof( ncls ); // histograms for profiles
468 
469  for (Int_t cls = 0; cls < ncls; cls++) {
470  hVars.at(cls).resize ( nvar+ntgt );
471  hVars.at(cls).assign ( nvar+ntgt, 0 ); // fill with zeros
472  mycorr.at(cls).resize( nvar+ntgt );
473  myprof.at(cls).resize( nvar+ntgt );
474  for (UInt_t ivar=0; ivar < nvar+ntgt; ivar++) {
475  mycorr.at(cls).at(ivar).resize( nvar+ntgt );
476  myprof.at(cls).at(ivar).resize( nvar+ntgt );
477  mycorr.at(cls).at(ivar).assign( nvar+ntgt, 0 ); // fill with zeros
478  myprof.at(cls).at(ivar).assign( nvar+ntgt, 0 ); // fill with zeros
479  }
480  }
481 
482  // if there are too many input variables, the creation of correlations plots blows up
483  // memory and basically kills the TMVA execution
484  // --> avoid above critical number (which can be user defined)
485  if (nvar+ntgt > (UInt_t)gConfig().GetVariablePlotting().fMaxNumOfAllowedVariablesForScatterPlots) {
486  Int_t nhists = (nvar+ntgt)*(nvar+ntgt - 1)/2;
487  Log() << kINFO << gTools().Color("dgreen") << Endl;
488  Log() << kINFO << "<PlotVariables> Will not produce scatter plots ==> " << Endl;
489  Log() << kINFO
490  << "| The number of " << nvar << " input variables and " << ntgt << " target values would require "
491  << nhists << " two-dimensional" << Endl;
492  Log() << kINFO
493  << "| histograms, which would occupy the computer's memory. Note that this" << Endl;
494  Log() << kINFO
495  << "| suppression does not have any consequences for your analysis, other" << Endl;
496  Log() << kINFO
497  << "| than not disposing of these scatter plots. You can modify the maximum" << Endl;
498  Log() << kINFO
499  << "| number of input variables allowed to generate scatter plots in your" << Endl;
500  Log() << "| script via the command line:" << Endl;
501  Log() << kINFO
502  << "| \"(TMVA::gConfig().GetVariablePlotting()).fMaxNumOfAllowedVariablesForScatterPlots = <some int>;\""
503  << gTools().Color("reset") << Endl;
504  Log() << Endl;
505  Log() << kINFO << "Some more output" << Endl;
506  }
507 
511 
512  for (UInt_t var_tgt = 0; var_tgt < 2; var_tgt++) { // create the histos first for the variables, then for the targets
513  UInt_t nloops = ( var_tgt == 0? nvar:ntgt ); // number of variables or number of targets
514  for (UInt_t ivar=0; ivar<nloops; ivar++) {
515  const VariableInfo& info = ( var_tgt == 0 ? Variable( ivar ) : Target(ivar) ); // choose the appropriate one (variable or target)
516  TString myVari = info.GetInternalName();
517 
518  Double_t mean = fVariableStats.at(fNumC-1).at( ( var_tgt*nvar )+ivar).fMean;
519  Double_t rms = fVariableStats.at(fNumC-1).at( ( var_tgt*nvar )+ivar).fRMS;
520 
521  for (Int_t cls = 0; cls < ncls; cls++) {
522 
523  TString className = fDataSetInfo.GetClassInfo(cls)->GetName();
524 
525  // add "target" in case of target variable (required for plotting macros)
526  className += (ntgt == 1 && var_tgt == 1 ? "_target" : "");
527 
528  // choose reasonable histogram ranges, by removing outliers
529  TH1* h = 0;
530  if (info.GetVarType() == 'I') {
531  // special treatment for integer variables
532  Int_t xmin = TMath::Nint( GetMin( ( var_tgt*nvar )+ivar) );
533  Int_t xmax = TMath::Nint( GetMax( ( var_tgt*nvar )+ivar) + 1 );
534  Int_t nbins = xmax - xmin;
535 
536  h = new TH1F( Form("%s__%s%s", myVari.Data(), className.Data(), transfType.Data()),
537  info.GetTitle(), nbins, xmin, xmax );
538  }
539  else {
540  Double_t xmin = TMath::Max( GetMin( ( var_tgt*nvar )+ivar), mean - timesRMS*rms );
541  Double_t xmax = TMath::Min( GetMax( ( var_tgt*nvar )+ivar), mean + timesRMS*rms );
542 
543  //std::cout << "Class="<<cls<<" xmin="<<xmin << " xmax="<<xmax<<" mean="<<mean<<" rms="<<rms<<" timesRMS="<<timesRMS<<std::endl;
544  // protection
545  if (xmin >= xmax) xmax = xmin*1.1; // try first...
546  if (xmin >= xmax) xmax = xmin + 1; // this if xmin == xmax == 0
547  // safety margin for values equal to the maximum within the histogram
548  xmax += (xmax - xmin)/nbins1D;
549 
550  h = new TH1F( Form("%s__%s%s", myVari.Data(), className.Data(), transfType.Data()),
551  info.GetTitle(), nbins1D, xmin, xmax );
552  }
553 
554  h->GetXaxis()->SetTitle( gTools().GetXTitleWithUnit( GetVariableAxisTitle( info ), info.GetUnit() ) );
555  h->GetYaxis()->SetTitle( gTools().GetYTitleWithUnit( *h, info.GetUnit(), kFALSE ) );
556  hVars.at(cls).at((var_tgt*nvar)+ivar) = h;
557 
558  // profile and scatter plots
559  if (nvar+ntgt <= (UInt_t)gConfig().GetVariablePlotting().fMaxNumOfAllowedVariablesForScatterPlots) {
560 
561  for (UInt_t v_t = 0; v_t < 2; v_t++) {
562  UInt_t nl = ( v_t==0?nvar:ntgt );
563  UInt_t start = ( v_t==0? (var_tgt==0?ivar+1:0):(var_tgt==0?nl:ivar+1) );
564  for (UInt_t j=start; j<nl; j++) {
565  // choose the appropriate one (variable or target)
566  const VariableInfo& infoj = ( v_t == 0 ? Variable( j ) : Target(j) );
567  TString myVarj = infoj.GetInternalName();
568 
569  Double_t rxmin = fVariableStats.at(fNumC-1).at( ( v_t*nvar )+ivar).fMin;
570  Double_t rxmax = fVariableStats.at(fNumC-1).at( ( v_t*nvar )+ivar).fMax;
571  Double_t rymin = fVariableStats.at(fNumC-1).at( ( v_t*nvar )+j).fMin;
572  Double_t rymax = fVariableStats.at(fNumC-1).at( ( v_t*nvar )+j).fMax;
573 
574  // scatter plot
575  TH2F* h2 = new TH2F( Form( "scat_%s_vs_%s_%s%s" , myVarj.Data(), myVari.Data(),
576  className.Data(), transfType.Data() ),
577  Form( "%s versus %s (%s)%s", infoj.GetTitle().Data(), info.GetTitle().Data(),
578  className.Data(), transfType.Data() ),
579  nbins2D, rxmin , rxmax,
580  nbins2D, rymin , rymax );
581 
582  h2->GetXaxis()->SetTitle( gTools().GetXTitleWithUnit( GetVariableAxisTitle( info ), info .GetUnit() ) );
583  h2->GetYaxis()->SetTitle( gTools().GetXTitleWithUnit( GetVariableAxisTitle( infoj ), infoj.GetUnit() ) );
584  mycorr.at(cls).at((var_tgt*nvar)+ivar).at((v_t*nvar)+j) = h2;
585 
586  // profile plot
587  TProfile* p = new TProfile( Form( "prof_%s_vs_%s_%s%s", myVarj.Data(),
588  myVari.Data(), className.Data(),
589  transfType.Data() ),
590  Form( "profile %s versus %s (%s)%s",
591  infoj.GetTitle().Data(), info.GetTitle().Data(),
592  className.Data(), transfType.Data() ), nbins1D,
593  rxmin, rxmax );
594  // info.GetMin(), info.GetMax() );
595 
596  p->GetXaxis()->SetTitle( gTools().GetXTitleWithUnit( GetVariableAxisTitle( info ), info .GetUnit() ) );
597  p->GetYaxis()->SetTitle( gTools().GetXTitleWithUnit( GetVariableAxisTitle( infoj ), infoj.GetUnit() ) );
598  myprof.at(cls).at((var_tgt*nvar)+ivar).at((v_t*nvar)+j) = p;
599  }
600  }
601  }
602  }
603  }
604  }
605 
606  UInt_t nevts = events.size();
607 
608  // compute correlation coefficient between target value and variables (regression only)
609  std::vector<Double_t> xregmean ( nvar+1, 0 );
610  std::vector<Double_t> x2regmean( nvar+1, 0 );
611  std::vector<Double_t> xCregmean( nvar+1, 0 );
612 
613  // fill the histograms (this approach should be faster than individual projection
614  for (UInt_t ievt=0; ievt<nevts; ievt++) {
615 
616  const Event* ev = events[ievt];
617 
618  Float_t weight = ev->GetWeight();
619  Int_t cls = ev->GetClass();
620 
621  // average correlation between first target and variables (so far only for single-target regression)
622  if (ntgt == 1) {
623  Float_t valr = ev->GetTarget(0);
624  xregmean[nvar] += valr;
625  x2regmean[nvar] += valr*valr;
626  for (UInt_t ivar=0; ivar<nvar; ivar++) {
627  Float_t vali = ev->GetValue(ivar);
628  xregmean[ivar] += vali;
629  x2regmean[ivar] += vali*vali;
630  xCregmean[ivar] += vali*valr;
631  }
632  }
633 
634  // fill correlation histograms
635  for (UInt_t var_tgt = 0; var_tgt < 2; var_tgt++) { // create the histos first for the variables, then for the targets
636  UInt_t nloops = ( var_tgt == 0? nvar:ntgt ); // number of variables or number of targets
637  for (UInt_t ivar=0; ivar<nloops; ivar++) {
638  Float_t vali = ( var_tgt == 0 ? ev->GetValue(ivar) : ev->GetTarget(ivar) );
639 
640  // variable histos
641  hVars.at(cls).at( ( var_tgt*nvar )+ivar)->Fill( vali, weight );
642 
643  // correlation histos
644  if (nvar+ntgt <= (UInt_t)gConfig().GetVariablePlotting().fMaxNumOfAllowedVariablesForScatterPlots) {
645 
646  for (UInt_t v_t = 0; v_t < 2; v_t++) {
647  UInt_t nl = ( v_t==0 ? nvar : ntgt );
648  UInt_t start = ( v_t==0 ? (var_tgt==0?ivar+1:0) : (var_tgt==0?nl:ivar+1) );
649  for (UInt_t j=start; j<nl; j++) {
650  Float_t valj = ( v_t == 0 ? ev->GetValue(j) : ev->GetTarget(j) );
651  mycorr.at(cls).at( ( var_tgt*nvar )+ivar).at( ( v_t*nvar )+j)->Fill( vali, valj, weight );
652  myprof.at(cls).at( ( var_tgt*nvar )+ivar).at( ( v_t*nvar )+j)->Fill( vali, valj, weight );
653  }
654  }
655  }
656  }
657  }
658  }
659 
660  // correlation analysis for ranking (single-target regression only)
661  if (ntgt == 1) {
662  for (UInt_t ivar=0; ivar<=nvar; ivar++) {
663  xregmean[ivar] /= nevts;
664  x2regmean[ivar] = x2regmean[ivar]/nevts - xregmean[ivar]*xregmean[ivar];
665  }
666  for (UInt_t ivar=0; ivar<nvar; ivar++) {
667  xCregmean[ivar] = xCregmean[ivar]/nevts - xregmean[ivar]*xregmean[nvar];
668  xCregmean[ivar] /= TMath::Sqrt( x2regmean[ivar]*x2regmean[nvar] );
669  }
670 
671  fRanking.push_back( new Ranking( GetName() + "Transformation", "|Correlation with target|" ) );
672  for (UInt_t ivar=0; ivar<nvar; ivar++) {
673  Double_t abscor = TMath::Abs( xCregmean[ivar] );
674  fRanking.back()->AddRank( Rank( fDataSetInfo.GetVariableInfo(ivar).GetLabel(), abscor ) );
675  }
676 
678 
679  // compute also mutual information (non-linear correlation measure)
680  fRanking.push_back( new Ranking( GetName() + "Transformation", "Mutual information" ) );
681  for (UInt_t ivar=0; ivar<nvar; ivar++) {
682  TH2F* h1 = mycorr.at(0).at( nvar ).at( ivar );
683  Double_t mi = gTools().GetMutualInformation( *h1 );
684  fRanking.back()->AddRank( Rank( fDataSetInfo.GetVariableInfo(ivar).GetLabel(), mi ) );
685  }
686 
687  // compute correlation ratio (functional correlations measure)
688  fRanking.push_back( new Ranking( GetName() + "Transformation", "Correlation Ratio" ) );
689  for (UInt_t ivar=0; ivar<nvar; ivar++) {
690  TH2F* h2 = mycorr.at(0).at( nvar ).at( ivar );
691  Double_t cr = gTools().GetCorrelationRatio( *h2 );
692  fRanking.back()->AddRank( Rank( fDataSetInfo.GetVariableInfo(ivar).GetLabel(), cr ) );
693  }
694 
695  // additionally compute correlation ratio from transposed histograms since correlation ratio is asymmetric
696  fRanking.push_back( new Ranking( GetName() + "Transformation", "Correlation Ratio (T)" ) );
697  for (UInt_t ivar=0; ivar<nvar; ivar++) {
698  TH2F* h2T = gTools().TransposeHist( *mycorr.at(0).at( nvar ).at( ivar ) );
699  Double_t cr = gTools().GetCorrelationRatio( *h2T );
700  fRanking.back()->AddRank( Rank( fDataSetInfo.GetVariableInfo(ivar).GetLabel(), cr ) );
701  delete h2T;
702  }
703  }
704  }
705  // computes ranking of input variables
706  // separation for 2-class classification
707  else if (fDataSetInfo.GetNClasses() == 2
708  && fDataSetInfo.GetClassInfo("Signal") != NULL
709  && fDataSetInfo.GetClassInfo("Background") != NULL
710  ) { // TODO: ugly hack.. adapt to new framework
711  fRanking.push_back( new Ranking( GetName() + "Transformation", "Separation" ) );
712  for (UInt_t i=0; i<nvar; i++) {
713  Double_t sep = gTools().GetSeparation( hVars.at(fDataSetInfo.GetClassInfo("Signal") ->GetNumber()).at(i),
714  hVars.at(fDataSetInfo.GetClassInfo("Background")->GetNumber()).at(i) );
715  fRanking.back()->AddRank( Rank( hVars.at(fDataSetInfo.GetClassInfo("Signal")->GetNumber()).at(i)->GetTitle(),
716  sep ) );
717  }
718  }
719 
720  // for regression compute performance from correlation with target value
721 
722  // write histograms
723 
724  TDirectory* localDir = theDirectory;
725  if (theDirectory == 0) {
726  // create directory in root dir
727  fRootBaseDir->cd();
728  TString outputDir = TString("InputVariables");
729  TListIter trIt(&fTransformations);
730  while (VariableTransformBase *trf = (VariableTransformBase*) trIt())
731  outputDir += "_" + TString(trf->GetShortName());
732 
733  TString uniqueOutputDir = outputDir;
734  Int_t counter = 0;
735  TObject* o = NULL;
736  while( (o = fRootBaseDir->FindObject(uniqueOutputDir)) != 0 ){
737  uniqueOutputDir = outputDir+Form("_%d",counter);
738  Log() << kINFO << "A " << o->ClassName() << " with name " << o->GetName() << " already exists in "
739  << fRootBaseDir->GetPath() << ", I will try with "<<uniqueOutputDir<<"." << Endl;
740  ++counter;
741  }
742 
743  // TObject* o = fRootBaseDir->FindObject(outputDir);
744  // if (o != 0) {
745  // Log() << kFATAL << "A " << o->ClassName() << " with name " << o->GetName() << " already exists in "
746  // << fRootBaseDir->GetPath() << "("<<outputDir<<")" << Endl;
747  // }
748  localDir = fRootBaseDir->mkdir( uniqueOutputDir );
749  localDir->cd();
750 
751  Log() << kVERBOSE << "Create and switch to directory " << localDir->GetPath() << Endl;
752  }
753  else {
754  theDirectory->cd();
755  }
756 
757  for (UInt_t i=0; i<nvar+ntgt; i++) {
758  for (Int_t cls = 0; cls < ncls; cls++) {
759  if (hVars.at(cls).at(i) != 0) {
760  hVars.at(cls).at(i)->Write();
761  hVars.at(cls).at(i)->SetDirectory(0);
762  delete hVars.at(cls).at(i);
763  }
764  }
765  }
766 
767  // correlation plots have dedicated directory
768  if (nvar+ntgt <= (UInt_t)gConfig().GetVariablePlotting().fMaxNumOfAllowedVariablesForScatterPlots) {
769 
770  localDir = localDir->mkdir( "CorrelationPlots" );
771  localDir ->cd();
772  Log() << kINFO << "Create scatter and profile plots in target-file directory: " << Endl;
773  Log() << kINFO << localDir->GetPath() << Endl;
774 
775 
776  for (UInt_t i=0; i<nvar+ntgt; i++) {
777  for (UInt_t j=i+1; j<nvar+ntgt; j++) {
778  for (Int_t cls = 0; cls < ncls; cls++) {
779  if (mycorr.at(cls).at(i).at(j) != 0 ) {
780  mycorr.at(cls).at(i).at(j)->Write();
781  mycorr.at(cls).at(i).at(j)->SetDirectory(0);
782  delete mycorr.at(cls).at(i).at(j);
783  }
784  if (myprof.at(cls).at(i).at(j) != 0) {
785  myprof.at(cls).at(i).at(j)->Write();
786  myprof.at(cls).at(i).at(j)->SetDirectory(0);
787  delete myprof.at(cls).at(i).at(j);
788  }
789  }
790  }
791  }
792  }
793  if (theDirectory != 0 ) theDirectory->cd();
794  else fRootBaseDir->cd();
795 }
796 
797 ////////////////////////////////////////////////////////////////////////////////
798 /// returns string for transformation
799 
801 {
802  VariableTransformBase* trf = ((VariableTransformBase*)GetTransformationList().Last());
803  if (!trf) return 0;
804  else return trf->GetTransformationStrings( fTransformationsReferenceClasses.back() );
805 }
806 
807 ////////////////////////////////////////////////////////////////////////////////
808 /// returns string for transformation
809 
811 {
812  VariableTransformBase* trf = ((VariableTransformBase*)GetTransformationList().Last());
813  if (!trf) return 0;
814  else return trf->GetName();
815 }
816 
817 ////////////////////////////////////////////////////////////////////////////////
818 /// write transformatino to stream
819 
820 void TMVA::TransformationHandler::WriteToStream( std::ostream& o ) const
821 {
822  TListIter trIt(&fTransformations);
823  std::vector< Int_t >::const_iterator rClsIt = fTransformationsReferenceClasses.begin();
824 
825  o << "NTransformtations " << fTransformations.GetSize() << std::endl << std::endl;
826 
827  ClassInfo* ci;
828  UInt_t i = 1;
829  while (VariableTransformBase *trf = (VariableTransformBase*) trIt()) {
830  o << "#TR -*-*-*-*-*-*-* transformation " << i++ << ": " << trf->GetName() << " -*-*-*-*-*-*-*-" << std::endl;
831  trf->WriteTransformationToStream(o);
832  ci = fDataSetInfo.GetClassInfo( (*rClsIt) );
833  TString clsName;
834  if (ci == 0 ) clsName = "AllClasses";
835  else clsName = ci->GetName();
836  o << "ReferenceClass " << clsName << std::endl;
837  rClsIt++;
838  }
839 }
840 
841 
842 ////////////////////////////////////////////////////////////////////////////////
843 /// XML node describing the transformation
844 /// return;
845 
846 void TMVA::TransformationHandler::AddXMLTo( void* parent ) const
847 {
848  if(!parent) return;
849  void* trfs = gTools().AddChild(parent, "Transformations");
850  gTools().AddAttr( trfs, "NTransformations", fTransformations.GetSize() );
851  TListIter trIt(&fTransformations);
852  while (VariableTransformBase *trf = (VariableTransformBase*) trIt()) trf->AttachXMLTo(trfs);
853 }
854 
855 ////////////////////////////////////////////////////////////////////////////////
856 ///VariableTransformBase* trf = ((VariableTransformBase*)GetTransformationList().Last());
857 ///trf->ReadTransformationFromStream(fin);
858 
860 {
861  Log() << kFATAL << "Read transformations not implemented" << Endl;
862  // TODO
863 }
864 
865 ////////////////////////////////////////////////////////////////////////////////
866 
868 {
869  void* ch = gTools().GetChild( trfsnode );
870  while(ch) {
871  Int_t idxCls = -1;
872  TString trfname;
873  gTools().ReadAttr(ch, "Name", trfname);
874 
875  VariableTransformBase* newtrf = 0;
876 
877  if (trfname == "Decorrelation" ) {
878  newtrf = new VariableDecorrTransform(fDataSetInfo);
879  }
880  else if (trfname == "PCA" ) {
881  newtrf = new VariablePCATransform(fDataSetInfo);
882  }
883  else if (trfname == "Gauss" ) {
884  newtrf = new VariableGaussTransform(fDataSetInfo);
885  }
886  else if (trfname == "Uniform" ) {
887  newtrf = new VariableGaussTransform(fDataSetInfo, "Uniform");
888  }
889  else if (trfname == "Normalize" ) {
890  newtrf = new VariableNormalizeTransform(fDataSetInfo);
891  }
892  else if (trfname == "Rearrange" ) {
893  newtrf = new VariableRearrangeTransform(fDataSetInfo);
894  }
895  else if (trfname != "None") {
896  }
897  else {
898  Log() << kFATAL << "<ReadFromXML> Variable transform '"
899  << trfname << "' unknown." << Endl;
900  }
901  newtrf->ReadFromXML( ch );
902  AddTransformation( newtrf, idxCls );
903  ch = gTools().GetNextChild(ch);
904  }
905 }
906 
907 ////////////////////////////////////////////////////////////////////////////////
908 /// prints ranking of input variables
909 
911 {
912  Log() << kINFO << " " << Endl;
913  Log() << kINFO << "Ranking input variables (method unspecific)..." << Endl;
914  std::vector<Ranking*>::const_iterator it = fRanking.begin();
915  for (; it != fRanking.end(); it++) (*it)->Print();
916 }
917 
918 ////////////////////////////////////////////////////////////////////////////////
919 
921 {
922  try {
923  return fVariableStats.at(cls).at(ivar).fMean;
924  }
925  catch(...) {
926  try {
927  return fVariableStats.at(fNumC-1).at(ivar).fMean;
928  }
929  catch(...) {
930  Log() << kWARNING << "Inconsistent variable state when reading the mean value. " << Endl;
931  }
932  }
933  Log() << kWARNING << "Inconsistent variable state when reading the mean value. Value 0 given back" << Endl;
934  return 0;
935 }
936 
937 
938 ////////////////////////////////////////////////////////////////////////////////
939 
941 {
942  try {
943  return fVariableStats.at(cls).at(ivar).fRMS;
944  }
945  catch(...) {
946  try {
947  return fVariableStats.at(fNumC-1).at(ivar).fRMS;
948  }
949  catch(...) {
950  Log() << kWARNING << "Inconsistent variable state when reading the RMS value. " << Endl;
951  }
952  }
953  Log() << kWARNING << "Inconsistent variable state when reading the RMS value. Value 0 given back" << Endl;
954  return 0;
955 }
956 
957 ////////////////////////////////////////////////////////////////////////////////
958 
960 {
961  try {
962  return fVariableStats.at(cls).at(ivar).fMin;
963  }
964  catch(...) {
965  try {
966  return fVariableStats.at(fNumC-1).at(ivar).fMin;
967  }
968  catch(...) {
969  Log() << kWARNING << "Inconsistent variable state when reading the minimum value. " << Endl;
970  }
971  }
972  Log() << kWARNING << "Inconsistent variable state when reading the minimum value. Value 0 given back" << Endl;
973  return 0;
974 }
975 
976 ////////////////////////////////////////////////////////////////////////////////
977 
979 {
980  try {
981  return fVariableStats.at(cls).at(ivar).fMax;
982  }
983  catch(...) {
984  try {
985  return fVariableStats.at(fNumC-1).at(ivar).fMax;
986  }
987  catch(...) {
988  Log() << kWARNING << "Inconsistent variable state when reading the maximum value. " << Endl;
989  }
990  }
991  Log() << kWARNING << "Inconsistent variable state when reading the maximum value. Value 0 given back" << Endl;
992  return 0;
993 }
std::vector< std::vector< TMVA::TransformationHandler::VariableStat > > fVariableStats
reference classes for the transformations
Int_t fMaxNumOfAllowedVariablesForScatterPlots
Definition: Config.h:88
float xmin
Definition: THbookFile.cxx:93
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
Ssiz_t Length() const
Definition: TString.h:390
float Float_t
Definition: RtypesCore.h:53
TString GetVariableAxisTitle(const VariableInfo &info) const
incorporates transformation type into title axis (usually for histograms)
TransformationHandler(DataSetInfo &, const TString &callerName)
constructor
UInt_t GetNClasses() const
Definition: DataSetInfo.h:152
Config & gConfig()
TH1 * h
Definition: legend2.C:5
UInt_t GetNTargets() const
Definition: DataSetInfo.h:129
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
Double_t GetMutualInformation(const TH2F &)
Mutual Information method for non-linear correlations estimates in 2D histogram Author: Moritz Backes...
Definition: Tools.cxx:598
void WriteToStream(std::ostream &o) const
write transformatino to stream
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
virtual TDirectory * mkdir(const char *name, const char *title="")
Create a sub-directory and return a pointer to the created directory.
Definition: TDirectory.cxx:955
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not...
Definition: Event.cxx:376
int nbins[3]
const std::vector< Event * > * CalcTransformations(const std::vector< Event * > &, Bool_t createNewVector=kFALSE)
computation of transformation
static std::string format(double x, double y, int digits, int width)
Profile Historam.
Definition: TProfile.h:34
void AddAttr(void *node, const char *, const T &value, Int_t precision=16)
Definition: Tools.h:308
UInt_t GetNVariables() const
Definition: DataSetInfo.h:128
virtual const char * GetPath() const
Returns the full path of the directory.
Definition: TDirectory.cxx:909
void * AddChild(void *parent, const char *childname, const char *content=0, bool isRootNode=false)
add child node
Definition: Tools.cxx:1134
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TObject Int_t at
Iterator of linked list.
Definition: TList.h:187
void MakeFunction(std::ostream &fout, const TString &fncName, Int_t part) const
create transformation function
const TString & GetInternalName() const
Definition: VariableInfo.h:63
void ReadFromStream(std::istream &istr)
VariableTransformBase* trf = ((VariableTransformBase*)GetTransformationList().Last()); trf->ReadTrans...
const char * Data() const
Definition: TString.h:349
Tools & gTools()
Definition: Tools.cxx:79
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
void SetTransformationReferenceClass(Int_t cls)
overrides the setting for all classes! (this is put in basically for the likelihood-method) be carefu...
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
Double_t GetRMS(Int_t ivar, Int_t cls=-1) const
void SetCallerName(const TString &name)
std::vector< std::vector< double > > Data
TH2D * h2
Definition: fit2dHist.C:45
TH1F * h1
Definition: legend1.C:5
virtual void ReadFromXML(void *trfnode)=0
Read the input variables from the XML node.
const Event * Transform(const Event *) const
the transformation
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:594
char GetVarType() const
Definition: VariableInfo.h:67
void PlotVariables(const std::vector< Event * > &events, TDirectory *theDirectory=0)
create histograms from the input variables
void AddStats(Int_t k, UInt_t ivar, Double_t mean, Double_t rms, Double_t min, Double_t max)
const TString & GetTitle() const
Definition: VariableInfo.h:65
const TString & GetUnit() const
Definition: VariableInfo.h:66
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:187
void CalcStats(const std::vector< Event * > &events)
virtual std::vector< TString > * GetTransformationStrings(Int_t cls) const
TODO –> adapt to variable,target,spectator selection default transformation output –> only indicate t...
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
Double_t GetSeparation(TH1 *S, TH1 *B) const
compute "separation" defined as <s2> = (1/2) Int_-oo..+oo { (S(x) - B(x))^2/(S(x) + B(x)) dx } ...
Definition: Tools.cxx:136
const TString & GetName() const
Definition: ClassInfo.h:72
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
TAxis * GetYaxis()
Definition: TH1.h:320
float xmax
Definition: THbookFile.cxx:93
void AddXMLTo(void *parent=0) const
XML node describing the transformation return;.
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:415
double Double_t
Definition: RtypesCore.h:55
Describe directory structure in memory.
Definition: TDirectory.h:44
virtual const char * GetName() const
Returns name of object.
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
Float_t GetValue(UInt_t ivar) const
return value of i'th variable
Definition: Event.cxx:231
The TH1 histogram class.
Definition: TH1.h:80
UInt_t GetClass() const
Definition: Event.h:86
std::vector< TString > * GetTransformationStringsOfLastTransform() const
returns string for transformation
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
#define name(a, b)
Definition: linkTestLib0.cpp:5
Double_t GetMin(Int_t ivar, Int_t cls=-1) const
void PrintVariableRanking() const
prints ranking of input variables
const TString & Color(const TString &)
human readable color strings
Definition: Tools.cxx:837
Mother of all ROOT objects.
Definition: TObject.h:58
Float_t GetTarget(UInt_t itgt) const
Definition: Event.h:101
Double_t GetMean(Int_t ivar, Int_t cls=-1) const
void SetSource(const std::string &source)
Definition: MsgLogger.h:74
virtual Bool_t cd(const char *path=0)
Change current directory to "this" directory.
Definition: TDirectory.cxx:433
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
#define NULL
Definition: Rtypes.h:82
TString()
TString default ctor.
Definition: TString.cxx:87
Double_t GetMax(Int_t ivar, Int_t cls=-1) const
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
TString GetName() const
return transformation name
const Bool_t kIterBackward
Definition: TCollection.h:44
TH2F * TransposeHist(const TH2F &)
Transpose quadratic histogram.
Definition: Tools.cxx:666
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:459
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition: TObject.cxx:379
Int_t Nint(T x)
Definition: TMath.h:480
VariablePlotting & GetVariablePlotting()
Definition: Config.h:77
const char * GetNameOfLastTransform() const
returns string for transformation
Double_t GetCorrelationRatio(const TH2F &)
Compute Correlation Ratio of 2D histogram to estimate functional dependency between two variables Aut...
Definition: Tools.cxx:629
Definition: math.cpp:60
TAxis * GetXaxis()
Definition: TH1.h:319
VariableTransformBase * AddTransformation(VariableTransformBase *, Int_t cls)
const Event * InverseTransform(const Event *, Bool_t suppressIfNoTargets=true) const