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