ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
VariableNormalizeTransform.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Joerg Stelzer, Helge Voss, Peter Speckmayer, Eckhard von Toerne
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : VariableNormalizeTransform *
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  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - CERN, Switzerland *
16  * Peter Speckmayer <Peter.Speckmayer@cern.ch> - CERN, Switzerland *
17  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
18  * Eckhard v. Toerne <evt@uni-bonn.de> - U of Bonn, Germany *
19  * *
20  * Copyright (c) 2005-2011: *
21  * CERN, Switzerland *
22  * MPI-K Heidelberg, Germany *
23  * U. of Bonn, Germany *
24  * *
25  * Redistribution and use in source and binary forms, with or without *
26  * modification, are permitted according to the terms listed in LICENSE *
27  * (http://tmva.sourceforge.net/LICENSE) *
28  **********************************************************************************/
29 
31 
32 #include "TMVA/DataSet.h"
33 #include "TMVA/DataSetInfo.h"
34 #include "TMVA/Event.h"
35 #include "TMVA/MsgLogger.h"
36 #include "TMVA/Tools.h"
37 #include "TMVA/Types.h"
38 #include "TMVA/VariableInfo.h"
39 
40 #include "TMatrixD.h"
41 #include "TMatrixDBase.h"
42 #include "TVectorD.h"
43 #include "TVectorF.h"
44 
45 #include <iostream>
46 #include <iomanip>
47 #include <cfloat>
48 
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 /// constructor
53 
54 TMVA::VariableNormalizeTransform::VariableNormalizeTransform( DataSetInfo& dsi )
55  : VariableTransformBase( dsi, Types::kNormalized, "Norm" )
56 {
57 }
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 
62 }
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 /// initialization of the normalization transformation
66 
68 {
69  UInt_t inputSize = fGet.size();
70  Int_t numC = GetNClasses()+1;
71  if (GetNClasses() <= 1 ) numC = 1;
72 
73  fMin.resize( numC );
74  fMax.resize( numC );
75  for (Int_t i=0; i<numC; i++) {
76  fMin.at(i).resize(inputSize);
77  fMax.at(i).resize(inputSize);
78  fMin.at(i).assign(inputSize, 0);
79  fMax.at(i).assign(inputSize, 0);
80  }
81 }
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// prepare transformation
85 
87 {
88  if (!IsEnabled() || IsCreated()) return kTRUE;
89 
90  Log() << kINFO << "Preparing the transformation." << Endl;
91 
92  Initialize();
93 
94  CalcNormalizationParams( events );
95 
96  SetCreated( kTRUE );
97 
98  return kTRUE;
99 }
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 
104 {
105  // apply the normalization transformation
106  if (!IsCreated()) Log() << kFATAL << "Transformation not yet created" << Endl;
107 
108  // if cls (the class chosen by the user) not existing,
109  // assume that he wants to have the matrix for all classes together.
110  // if (cls < 0 || cls > GetNClasses()) {
111  // if (GetNClasses() > 1 ) cls = GetNClasses();
112  // else cls = (fMin.size()==1?0:2);
113  // }
114  // EVT this is a workaround to address the reader problem with transforma and EvaluateMVA(std::vector<float/double> ,...)
115  if (cls < 0 || cls >= (int) fMin.size()) cls = fMin.size()-1;
116  // EVT workaround end
117 
118  FloatVector input; // will be filled with the selected variables, targets, (spectators)
119  FloatVector output; // will be filled with the selected variables, targets, (spectators)
120  std::vector<Char_t> mask; // entries with kTRUE must not be transformed
121  GetInput( ev, input, mask );
122 
123  if (fTransformedEvent==0) fTransformedEvent = new Event();
124 
125  Float_t min,max;
126  const FloatVector& minVector = fMin.at(cls);
127  const FloatVector& maxVector = fMax.at(cls);
128 
129  UInt_t iidx = 0;
130  std::vector<Char_t>::iterator itMask = mask.begin();
131  for ( std::vector<Float_t>::iterator itInp = input.begin(), itInpEnd = input.end(); itInp != itInpEnd; ++itInp) { // loop over input variables
132  if( (*itMask) ){
133  ++iidx;
134  ++itMask;
135  // don't put any value into output if the value is masked
136  continue;
137  }
138 
139  Float_t val = (*itInp);
140 
141  min = minVector.at(iidx);
142  max = maxVector.at(iidx);
143  Float_t offset = min;
144  Float_t scale = 1.0/(max-min);
145 
146  Float_t valnorm = (val-offset)*scale * 2 - 1;
147  output.push_back( valnorm );
148 
149  ++iidx;
150  ++itMask;
151  }
152 
153  SetOutput( fTransformedEvent, output, mask, ev );
154  return fTransformedEvent;
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// apply the inverse transformation
159 
161 {
162  if (!IsCreated()) Log() << kFATAL << "Transformation not yet created" << Endl;
163 
164  // if cls (the class chosen by the user) not existing,
165  // assume that user wants to have the transformation for all classes together.
166  if (cls < 0 || cls > GetNClasses()) {
167  if (GetNClasses() > 1 ) cls = GetNClasses();
168  else cls = 0;
169  }
170 
171  FloatVector input; // will be filled with the selected variables, targets, (spectators)
172  FloatVector output; // will be filled with the output
173  std::vector<Char_t> mask;
174  GetInput( ev, input, mask, kTRUE );
175 
176  if (fBackTransformedEvent==0) fBackTransformedEvent = new Event( *ev );
177 
178  Float_t min,max;
179  const FloatVector& minVector = fMin.at(cls);
180  const FloatVector& maxVector = fMax.at(cls);
181 
182  UInt_t iidx = 0;
183  for ( std::vector<Float_t>::iterator itInp = input.begin(), itInpEnd = input.end(); itInp != itInpEnd; ++itInp) { // loop over input variables
184  Float_t val = (*itInp);
185 
186  min = minVector.at(iidx);
187  max = maxVector.at(iidx);
188  Float_t offset = min;
189  Float_t scale = 1.0/(max-min);
190 
191  Float_t valnorm = offset+((val+1)/(scale * 2));
192  output.push_back( valnorm );
193 
194  ++iidx;
195  }
196 
197  SetOutput( fBackTransformedEvent, output, mask, ev, kTRUE );
198 
199  return fBackTransformedEvent;
200 }
201 
202 ////////////////////////////////////////////////////////////////////////////////
203 /// compute offset and scale from min and max
204 
205 void TMVA::VariableNormalizeTransform::CalcNormalizationParams( const std::vector< Event*>& events )
206 {
207  if (events.size() <= 1)
208  Log() << kFATAL << "Not enough events (found " << events.size() << ") to calculate the normalization" << Endl;
209 
210  FloatVector input; // will be filled with the selected variables, targets, (spectators)
211  std::vector<Char_t> mask;
212 
213  UInt_t inputSize = fGet.size(); // number of input variables
214 
215  const UInt_t nCls = GetNClasses();
216  Int_t numC = nCls+1; // prepare the min and max values for each of the classes and additionally for all classes (if more than one)
217  Int_t all = nCls; // at idx the min and max values for "all" classes are stored
218  if (nCls <= 1 ) {
219  numC = 1;
220  all = 0;
221  }
222 
223  for (UInt_t iinp=0; iinp<inputSize; ++iinp) {
224  for (Int_t ic = 0; ic < numC; ic++) {
225  fMin.at(ic).at(iinp) = FLT_MAX;
226  fMax.at(ic).at(iinp) = -FLT_MAX;
227  }
228  }
229 
230  std::vector<Event*>::const_iterator evIt = events.begin();
231  for (;evIt!=events.end();evIt++) { // loop over all events
232  const TMVA::Event* event = (*evIt); // get the event
233 
234  UInt_t cls = (*evIt)->GetClass(); // get the class of this event
235 
236  FloatVector& minVector = fMin.at(cls);
237  FloatVector& maxVector = fMax.at(cls);
238 
239  FloatVector& minVectorAll = fMin.at(all);
240  FloatVector& maxVectorAll = fMax.at(all);
241 
242  GetInput(event,input,mask); // select the input variables for the transformation and get them from the event
243  UInt_t iidx = 0;
244  for ( std::vector<Float_t>::iterator itInp = input.begin(), itInpEnd = input.end(); itInp != itInpEnd; ++itInp) { // loop over input variables
245  Float_t val = (*itInp);
246 
247  if( minVector.at(iidx) > val ) minVector.at(iidx) = val;
248  if( maxVector.at(iidx) < val ) maxVector.at(iidx) = val;
249 
250  if (nCls != 1) { // in case more than one class exists, compute min and max as well for all classes together
251  if (minVectorAll.at(iidx) > val) minVectorAll.at(iidx) = val;
252  if (maxVectorAll.at(iidx) < val) maxVectorAll.at(iidx) = val;
253  }
254 
255  ++iidx;
256  }
257  }
258 
259  return;
260 }
261 
262 ////////////////////////////////////////////////////////////////////////////////
263 /// creates string with variable transformations applied
264 
266 {
267  // if cls (the class chosen by the user) not existing, assume that user wants to
268  // have the matrix for all classes together.
269  if (cls < 0 || cls > GetNClasses()) cls = GetNClasses();
270 
271  Float_t min, max;
272  const UInt_t size = fGet.size();
273  std::vector<TString>* strVec = new std::vector<TString>(size);
274 
275  UInt_t iinp = 0;
276  for( ItVarTypeIdxConst itGet = fGet.begin(), itGetEnd = fGet.end(); itGet != itGetEnd; ++itGet ) {
277  min = fMin.at(cls).at(iinp);
278  max = fMax.at(cls).at(iinp);
279 
280  Char_t type = (*itGet).first;
281  UInt_t idx = (*itGet).second;
282  Float_t offset = min;
283  Float_t scale = 1.0/(max-min);
284  TString str("");
285  VariableInfo& varInfo = (type=='v'?fDsi.GetVariableInfo(idx):(type=='t'?fDsi.GetTargetInfo(idx):fDsi.GetSpectatorInfo(idx)));
286 
287  if (offset < 0) str = Form( "2*%g*([%s] + %g) - 1", scale, varInfo.GetLabel().Data(), -offset );
288  else str = Form( "2*%g*([%s] - %g) - 1", scale, varInfo.GetLabel().Data(), offset );
289  (*strVec)[iinp] = str;
290 
291  ++iinp;
292  }
293 
294  return strVec;
295 }
296 
297 ////////////////////////////////////////////////////////////////////////////////
298 /// write the transformation to the stream
299 
301 {
302  o << "# min max for all variables for all classes one after the other and as a last entry for all classes together" << std::endl;
303 
304  Int_t numC = GetNClasses()+1;
305  if (GetNClasses() <= 1 ) numC = 1;
306 
307  UInt_t nvars = GetNVariables();
308  UInt_t ntgts = GetNTargets();
309 
310  for (Int_t icls = 0; icls < numC; icls++ ) {
311  o << icls << std::endl;
312  for (UInt_t ivar=0; ivar<nvars; ivar++)
313  o << std::setprecision(12) << std::setw(20) << fMin.at(icls).at(ivar) << " "
314  << std::setprecision(12) << std::setw(20) << fMax.at(icls).at(ivar) << std::endl;
315  for (UInt_t itgt=0; itgt<ntgts; itgt++)
316  o << std::setprecision(12) << std::setw(20) << fMin.at(icls).at(nvars+itgt) << " "
317  << std::setprecision(12) << std::setw(20) << fMax.at(icls).at(nvars+itgt) << std::endl;
318  }
319  o << "##" << std::endl;
320 }
321 
322 ////////////////////////////////////////////////////////////////////////////////
323 /// create XML description of Normalize transformation
324 
326 {
327  void* trfxml = gTools().AddChild(parent, "Transform");
328  gTools().AddAttr(trfxml, "Name", "Normalize");
330 
331  Int_t numC = (GetNClasses()<= 1)?1:GetNClasses()+1;
332 
333  for( Int_t icls=0; icls<numC; icls++ ) {
334  void* clsxml = gTools().AddChild(trfxml, "Class");
335  gTools().AddAttr(clsxml, "ClassIndex", icls);
336  void* inpxml = gTools().AddChild(clsxml, "Ranges");
337  UInt_t iinp = 0;
338  for( ItVarTypeIdx itGet = fGet.begin(), itGetEnd = fGet.end(); itGet != itGetEnd; ++itGet ) {
339  void* mmxml = gTools().AddChild(inpxml, "Range");
340  gTools().AddAttr(mmxml, "Index", iinp);
341  gTools().AddAttr(mmxml, "Min", fMin.at(icls).at(iinp) );
342  gTools().AddAttr(mmxml, "Max", fMax.at(icls).at(iinp) );
343  ++iinp;
344  }
345  }
346 }
347 
348 ////////////////////////////////////////////////////////////////////////////////
349 /// Read the transformation matrices from the xml node
350 
352 {
353  Bool_t newFormat = kFALSE;
354 
355  void* inpnode = NULL;
356 
357  inpnode = gTools().GetChild(trfnode, "Selection"); // new xml format
358  if( inpnode != NULL )
359  newFormat = kTRUE;
360 
361  if( newFormat ){
362  // ------------- new format --------------------
363  // read input
365 
366  // read transformation information
367 
368  UInt_t size = fGet.size();
369  UInt_t classindex, idx;
370 
371  void* ch = gTools().GetChild( trfnode, "Class" );
372  while(ch) {
373  Int_t ci = 0;
374  gTools().ReadAttr(ch, "ClassIndex", ci);
375  classindex = UInt_t(ci);
376 
377  fMin.resize(classindex+1);
378  fMax.resize(classindex+1);
379 
380  fMin[classindex].resize(size,Float_t(0));
381  fMax[classindex].resize(size,Float_t(0));
382 
383  void* clch = gTools().GetChild( ch );
384  while(clch) {
385  TString nodeName(gTools().GetName(clch));
386  if(nodeName=="Ranges") {
387  void* varch = gTools().GetChild( clch );
388  while(varch) {
389  gTools().ReadAttr(varch, "Index", idx);
390  gTools().ReadAttr(varch, "Min", fMin[classindex][idx]);
391  gTools().ReadAttr(varch, "Max", fMax[classindex][idx]);
392  varch = gTools().GetNextChild( varch );
393  }
394  }
395  clch = gTools().GetNextChild( clch );
396  }
397  ch = gTools().GetNextChild( ch );
398  }
399  SetCreated();
400  return;
401  }
402 
403  // ------------- old format --------------------
404  UInt_t classindex, varindex, tgtindex, nvars, ntgts;
405  // coverity[tainted_data_argument]
406  gTools().ReadAttr(trfnode, "NVariables", nvars);
407  // coverity[tainted_data_argument]
408  gTools().ReadAttr(trfnode, "NTargets", ntgts);
409  // coverity[tainted_data_argument]
410 
411  for( UInt_t ivar = 0; ivar < nvars; ++ivar ){
412  fGet.push_back(std::pair<Char_t,UInt_t>('v',ivar));
413  }
414  for( UInt_t itgt = 0; itgt < ntgts; ++itgt ){
415  fGet.push_back(std::pair<Char_t,UInt_t>('t',itgt));
416  }
417  void* ch = gTools().GetChild( trfnode );
418  while(ch) {
419  gTools().ReadAttr(ch, "ClassIndex", classindex);
420 
421  fMin.resize(classindex+1);
422  fMax.resize(classindex+1);
423  fMin[classindex].resize(nvars+ntgts,Float_t(0));
424  fMax[classindex].resize(nvars+ntgts,Float_t(0));
425 
426  void* clch = gTools().GetChild( ch );
427  while(clch) {
428  TString nodeName(gTools().GetName(clch));
429  if(nodeName=="Variables") {
430  void* varch = gTools().GetChild( clch );
431  while(varch) {
432  gTools().ReadAttr(varch, "VarIndex", varindex);
433  gTools().ReadAttr(varch, "Min", fMin[classindex][varindex]);
434  gTools().ReadAttr(varch, "Max", fMax[classindex][varindex]);
435  varch = gTools().GetNextChild( varch );
436  }
437  } else if (nodeName=="Targets") {
438  void* tgtch = gTools().GetChild( clch );
439  while(tgtch) {
440  gTools().ReadAttr(tgtch, "TargetIndex", tgtindex);
441  gTools().ReadAttr(tgtch, "Min", fMin[classindex][nvars+tgtindex]);
442  gTools().ReadAttr(tgtch, "Max", fMax[classindex][nvars+tgtindex]);
443  tgtch = gTools().GetNextChild( tgtch );
444  }
445  }
446  clch = gTools().GetNextChild( clch );
447  }
448  ch = gTools().GetNextChild( ch );
449  }
450  SetCreated();
451 }
452 
453 ////////////////////////////////////////////////////////////////////////////////
454 /// this method is only used when building a normalization transformation
455 /// from old text files
456 /// in this case regression didn't exist and there were no targets
457 
458 void TMVA::VariableNormalizeTransform::BuildTransformationFromVarInfo( const std::vector<TMVA::VariableInfo>& var )
459 {
460  UInt_t nvars = GetNVariables();
461 
462  if(var.size() != nvars)
463  Log() << kFATAL << "<BuildTransformationFromVarInfo> can't build transformation,"
464  << " since the number of variables disagree" << Endl;
465 
466  UInt_t numC = (GetNClasses()<=1)?1:GetNClasses()+1;
467  fMin.clear();fMin.resize( numC );
468  fMax.clear();fMax.resize( numC );
469 
470 
471  for(UInt_t cls=0; cls<numC; ++cls) {
472  fMin[cls].resize(nvars+GetNTargets(),0);
473  fMax[cls].resize(nvars+GetNTargets(),0);
474  UInt_t vidx(0);
475  for(std::vector<TMVA::VariableInfo>::const_iterator v = var.begin(); v!=var.end(); ++v, ++vidx) {
476  fMin[cls][vidx] = v->GetMin();
477  fMax[cls][vidx] = v->GetMax();
478  fGet.push_back(std::pair<Char_t,UInt_t>('v',vidx));
479  }
480  }
481  SetCreated();
482 }
483 
484 ////////////////////////////////////////////////////////////////////////////////
485 /// Read the variable ranges from an input stream
486 
488 {
489  UInt_t nvars = GetNVariables();
490  UInt_t ntgts = GetNTargets();
491  for( UInt_t ivar = 0; ivar < nvars; ++ivar ){
492  fGet.push_back(std::pair<Char_t,UInt_t>('v',ivar));
493  }
494  for( UInt_t itgt = 0; itgt < ntgts; ++itgt ){
495  fGet.push_back(std::pair<Char_t,UInt_t>('t',itgt));
496  }
497  char buf[512];
498  char buf2[512];
499  istr.getline(buf,512);
500  TString strvar, dummy;
501  Int_t icls;
502  TString test;
503  while (!(buf[0]=='#'&& buf[1]=='#')) { // if line starts with ## return
504  char* p = buf;
505  while (*p==' ' || *p=='\t') p++; // 'remove' leading whitespace
506  if (*p=='#' || *p=='\0') {
507  istr.getline(buf,512);
508  continue; // if comment or empty line, read the next line
509  }
510  std::stringstream sstr(buf);
511  sstr >> icls;
512  for (UInt_t ivar=0;ivar<nvars;ivar++) {
513  istr.getline(buf2,512); // reading the next line
514  std::stringstream sstr2(buf2);
515  sstr2 >> fMin[icls][ivar] >> fMax[icls][ivar];
516  }
517  for (UInt_t itgt=0;itgt<ntgts;itgt++) {
518  istr.getline(buf2,512); // reading the next line
519  std::stringstream sstr2(buf2);
520  sstr2 >> fMin[icls][nvars+itgt] >> fMax[icls][nvars+itgt];
521  }
522  istr.getline(buf,512); // reading the next line
523  }
524  SetCreated();
525 }
526 
527 ////////////////////////////////////////////////////////////////////////////////
528 /// prints the transformation ranges
529 
531 {
532  Int_t nCls = GetNClasses();
533  Int_t numC = nCls+1;
534  if (nCls <= 1 ) numC = 1;
535  for (Int_t icls = 0; icls < numC; icls++ ) {
536  if( icls == nCls )
537  Log() << kINFO << "Transformation for all classes based on these ranges:" << Endl;
538  else
539  Log() << kINFO << "Transformation for class " << icls << " based on these ranges:" << Endl;
540  UInt_t iinp = 0;
541  for( ItVarTypeIdxConst itGet = fGet.begin(), itGetEnd = fGet.end(); itGet != itGetEnd; ++itGet ){
542  Char_t type = (*itGet).first;
543  UInt_t idx = (*itGet).second;
544 
545  TString typeString = (type=='v'?"Variable: ": (type=='t'?"Target : ":"Spectator : ") );
546  Log() << typeString.Data() << std::setw(20) << fMin[icls][idx] << std::setw(20) << fMax[icls][idx] << Endl;
547 
548  ++iinp;
549  }
550  }
551 }
552 
553 ////////////////////////////////////////////////////////////////////////////////
554 /// creates a normalizing function
555 /// TODO include target-transformation into makefunction
556 
557 void TMVA::VariableNormalizeTransform::MakeFunction( std::ostream& fout, const TString& fcncName,
558  Int_t part, UInt_t trCounter, Int_t )
559 {
560  UInt_t nVar = fGet.size();
561  UInt_t numC = fMin.size();
562  if (part==1) {
563  fout << std::endl;
564  fout << " double fMin_"<<trCounter<<"["<<numC<<"]["<<nVar<<"];" << std::endl;
565  fout << " double fMax_"<<trCounter<<"["<<numC<<"]["<<nVar<<"];" << std::endl;
566  }
567 
568  if (part==2) {
569  fout << std::endl;
570  fout << "//_______________________________________________________________________" << std::endl;
571  fout << "inline void " << fcncName << "::InitTransform_"<<trCounter<<"()" << std::endl;
572  fout << "{" << std::endl;
573  fout << " // Normalization transformation, initialisation" << std::endl;
574  for (UInt_t ivar=0; ivar<nVar; ivar++) {
575  for (UInt_t icls = 0; icls < numC; icls++) {
576  Double_t min = TMath::Min( FLT_MAX, fMin.at(icls).at(ivar) );
577  Double_t max = TMath::Max(-FLT_MAX, fMax.at(icls).at(ivar) );
578  fout << " fMin_"<<trCounter<<"["<<icls<<"]["<<ivar<<"] = " << std::setprecision(12)
579  << min << ";" << std::endl;
580  fout << " fMax_"<<trCounter<<"["<<icls<<"]["<<ivar<<"] = " << std::setprecision(12)
581  << max << ";" << std::endl;
582  }
583  }
584  fout << "}" << std::endl;
585  fout << std::endl;
586  fout << "//_______________________________________________________________________" << std::endl;
587  fout << "inline void " << fcncName << "::Transform_"<<trCounter<<"( std::vector<double>& iv, int cls) const" << std::endl;
588  fout << "{" << std::endl;
589  fout << " // Normalization transformation" << std::endl;
590  fout << " if (cls < 0 || cls > "<<GetNClasses()<<") {"<< std::endl;
591  fout << " if ("<<GetNClasses()<<" > 1 ) cls = "<<GetNClasses()<<";"<< std::endl;
592  fout << " else cls = "<<(fMin.size()==1?0:2)<<";"<< std::endl;
593  fout << " }"<< std::endl;
594  fout << " const int nVar = " << nVar << ";" << std::endl << std::endl;
595  fout << " // get indices of used variables" << std::endl;
596  VariableTransformBase::MakeFunction(fout, fcncName, 0, trCounter, 0 );
597  fout << " static std::vector<double> dv;" << std::endl; // simply made it static so it doesn't need to be re-booked every time
598  fout << " dv.resize(nVar);" << std::endl;
599  fout << " for (int ivar=0; ivar<nVar; ivar++) dv[ivar] = iv[indicesGet.at(ivar)];" << std::endl;
600 
601  fout << " for (int ivar=0;ivar<"<<nVar<<";ivar++) {" << std::endl;
602  fout << " double offset = fMin_"<<trCounter<<"[cls][ivar];" << std::endl;
603  fout << " double scale = 1.0/(fMax_"<<trCounter<<"[cls][ivar]-fMin_"<<trCounter<<"[cls][ivar]);" << std::endl;
604  fout << " iv[indicesPut.at(ivar)] = (dv[ivar]-offset)*scale * 2 - 1;" << std::endl;
605  fout << " }" << std::endl;
606  fout << "}" << std::endl;
607  }
608 }
void CalcNormalizationParams(const std::vector< Event * > &events)
compute offset and scale from min and max
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
virtual void MakeFunction(std::ostream &fout, const TString &fncName, Int_t part, UInt_t trCounter, Int_t cls)
creates a normalizing function TODO include target-transformation into makefunction ...
float Float_t
Definition: RtypesCore.h:53
tuple offset
Definition: tree.py:93
virtual void MakeFunction(std::ostream &fout, const TString &fncName, Int_t part, UInt_t trCounter, Int_t cls)=0
getinput and setoutput equivalent
virtual void AttachXMLTo(void *parent)=0
create XML description the transformation (write out info of selected variables)
Basic string class.
Definition: TString.h:137
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Bool_t PrepareTransformation(const std::vector< Event * > &)
prepare transformation
void AddAttr(void *node, const char *, const T &value, Int_t precision=16)
Definition: Tools.h:308
void * AddChild(void *parent, const char *childname, const char *content=0, bool isRootNode=false)
add child node
Definition: Tools.cxx:1134
const char * Data() const
Definition: TString.h:349
Tools & gTools()
Definition: Tools.cxx:79
ClassImp(TMVA::VariableNormalizeTransform) TMVA
constructor
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
virtual const Event * Transform(const Event *const, Int_t cls) const
VectorOfCharAndInt::iterator ItVarTypeIdx
virtual void PrintTransformation(std::ostream &o)
prints the transformation ranges
virtual void ReadFromXML(void *trfnode)=0
Read the input variables from the XML node.
virtual const Event * InverseTransform(const Event *const, Int_t cls) const
apply the inverse transformation
void WriteTransformationToStream(std::ostream &) const
write the transformation to the stream
VectorOfCharAndInt::const_iterator ItVarTypeIdxConst
virtual void AttachXMLTo(void *parent)
create XML description of Normalize transformation
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
SVector< double, 2 > v
Definition: Dict.h:5
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
std::vector< TString > * GetTransformationStrings(Int_t cls) const
creates string with variable transformations applied
void Initialize()
initialization of the normalization transformation
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
void BuildTransformationFromVarInfo(const std::vector< TMVA::VariableInfo > &var)
this method is only used when building a normalization transformation from old text files in this cas...
double Double_t
Definition: RtypesCore.h:55
int type
Definition: TGX11.cxx:120
static RooMathCoreReg dummy
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
virtual void ReadFromXML(void *trfnode)
Read the transformation matrices from the xml node.
UInt_t GetClass() const
Definition: Event.h:86
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
char Char_t
Definition: RtypesCore.h:29
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
#define NULL
Definition: Rtypes.h:82
const TString & GetLabel() const
Definition: VariableInfo.h:64
static void output(int code)
Definition: gifencode.c:226
const Bool_t kTRUE
Definition: Rtypes.h:91
tuple all
Definition: na49view.py:13
Definition: math.cpp:60
void ReadTransformationFromStream(std::istream &, const TString &)
Read the variable ranges from an input stream.