ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TMLPAnalyzer.cxx
Go to the documentation of this file.
1 // @(#)root/mlp:$Id$
2 // Author: Christophe.Delaere@cern.ch 25/04/04
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2003, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 ///////////////////////////////////////////////////////////////////////////
13 //
14 // TMLPAnalyzer
15 //
16 // This utility class contains a set of tests usefull when developing
17 // a neural network.
18 // It allows you to check for unneeded variables, and to control
19 // the network structure.
20 //
21 ///////////////////////////////////////////////////////////////////////////
22 
23 #include "TROOT.h"
24 #include "TSynapse.h"
25 #include "TNeuron.h"
26 #include "TMultiLayerPerceptron.h"
27 #include "TMLPAnalyzer.h"
28 #include "TTree.h"
29 #include "TTreeFormula.h"
30 #include "TEventList.h"
31 #include "TH1D.h"
32 #include "TProfile.h"
33 #include "THStack.h"
34 #include "TLegend.h"
35 #include "TPad.h"
36 #include "TCanvas.h"
37 #include "TGaxis.h"
38 #include "TRegexp.h"
39 #include "TMath.h"
40 #include "Riostream.h"
41 #include <stdlib.h>
42 
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// Destructor
47 
49 {
50  delete fAnalysisTree;
51  delete fIOTree;
52 }
53 
54 ////////////////////////////////////////////////////////////////////////////////
55 /// Returns the number of layers.
56 
58 {
59  TString fStructure = fNetwork->GetStructure();
60  return fStructure.CountChar(':')+1;
61 }
62 
63 ////////////////////////////////////////////////////////////////////////////////
64 /// Returns the number of neurons in given layer.
65 
67 {
68  if(layer==1) {
69  TString fStructure = fNetwork->GetStructure();
70  TString input = TString(fStructure(0, fStructure.First(':')));
71  return input.CountChar(',')+1;
72  }
73  else if(layer==GetLayers()) {
74  TString fStructure = fNetwork->GetStructure();
75  TString output = TString(fStructure(fStructure.Last(':') + 1,
76  fStructure.Length() - fStructure.Last(':')));
77  return output.CountChar(',')+1;
78  }
79  else {
80  Int_t cnt=1;
81  TString fStructure = fNetwork->GetStructure();
82  TString hidden = TString(fStructure(fStructure.First(':') + 1,
83  fStructure.Last(':') - fStructure.First(':') - 1));
84  Int_t beg = 0;
85  Int_t end = hidden.Index(":", beg + 1);
86  Int_t num = 0;
87  while (end != -1) {
88  num = atoi(TString(hidden(beg, end - beg)).Data());
89  cnt++;
90  beg = end + 1;
91  end = hidden.Index(":", beg + 1);
92  if(layer==cnt) return num;
93  }
94  num = atoi(TString(hidden(beg, hidden.Length() - beg)).Data());
95  cnt++;
96  if(layer==cnt) return num;
97  }
98  return -1;
99 }
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 /// Returns the formula used as input for neuron (idx) in
103 /// the first layer.
104 
106 {
107  TString fStructure = fNetwork->GetStructure();
108  TString input = TString(fStructure(0, fStructure.First(':')));
109  Int_t beg = 0;
110  Int_t end = input.Index(",", beg + 1);
111  TString brName;
112  Int_t cnt = 0;
113  while (end != -1) {
114  brName = TString(input(beg, end - beg));
115  if (brName[0]=='@')
116  brName = brName(1,brName.Length()-1);
117  beg = end + 1;
118  end = input.Index(",", beg + 1);
119  if(cnt==idx) return brName;
120  cnt++;
121  }
122  brName = TString(input(beg, input.Length() - beg));
123  if (brName[0]=='@')
124  brName = brName(1,brName.Length()-1);
125  return brName;
126 }
127 
128 ////////////////////////////////////////////////////////////////////////////////
129 /// Returns the name of any neuron from the input layer
130 
132 {
133  TNeuron* neuron=(TNeuron*)fNetwork->fFirstLayer[in];
134  return neuron ? neuron->GetName() : "NO SUCH NEURON";
135 }
136 
137 ////////////////////////////////////////////////////////////////////////////////
138 /// Returns the name of any neuron from the output layer
139 
141 {
142  TNeuron* neuron=(TNeuron*)fNetwork->fLastLayer[out];
143  return neuron ? neuron->GetName() : "NO SUCH NEURON";
144 }
145 
146 ////////////////////////////////////////////////////////////////////////////////
147 /// Gives some information about the network in the terminal.
148 
150 {
151  TString fStructure = fNetwork->GetStructure();
152  std::cout << "Network with structure: " << fStructure.Data() << std::endl;
153  std::cout << "inputs with low values in the differences plot may not be needed" << std::endl;
154  // Checks if some input variable is not needed
155  char var[64], sel[64];
156  for (Int_t i = 0; i < GetNeurons(1); i++) {
157  snprintf(var,64,"diff>>tmp%d",i);
158  snprintf(sel,64,"inNeuron==%d",i);
159  fAnalysisTree->Draw(var, sel, "goff");
160  TH1F* tmp = (TH1F*)gDirectory->Get(Form("tmp%d",i));
161  if (!tmp) continue;
162  std::cout << GetInputNeuronTitle(i)
163  << " -> " << tmp->GetMean()
164  << " +/- " << tmp->GetRMS() << std::endl;
165  }
166 }
167 
168 ////////////////////////////////////////////////////////////////////////////////
169 /// Collect information about what is usefull in the network.
170 /// This method has to be called first when analyzing a network.
171 /// Fills the two analysis trees.
172 
174 {
175  Double_t shift = 0.1;
176  TTree* data = fNetwork->fData;
177  TEventList* test = fNetwork->fTest;
178  Int_t nEvents = test->GetN();
179  Int_t nn = GetNeurons(1);
180  Double_t* params = new Double_t[nn];
181  Double_t* rms = new Double_t[nn];
182  TTreeFormula** formulas = new TTreeFormula*[nn];
183  Int_t* index = new Int_t[nn];
184  TString formula;
185  TRegexp re("{[0-9]+}$");
186  Ssiz_t len = formula.Length();
187  Ssiz_t pos = -1;
188  Int_t i(0), j(0), k(0), l(0);
189  for(i=0; i<nn; i++){
190  formula = GetNeuronFormula(i);
191  pos = re.Index(formula,&len);
192  if(pos==-1 || len<3) {
193  formulas[i] = new TTreeFormula(Form("NF%lu",(ULong_t)this),formula,data);
194  index[i] = 0;
195  }
196  else {
197  TString newformula(formula,pos);
198  TString val = formula(pos+1,len-2);
199  formulas[i] = new TTreeFormula(Form("NF%lu",(ULong_t)this),newformula,data);
200  formula = newformula;
201  index[i] = val.Atoi();
202  }
203  TH1D tmp("tmpb", "tmpb", 1, -FLT_MAX, FLT_MAX);
204  data->Draw(Form("%s>>tmpb",formula.Data()),"","goff");
205  rms[i] = tmp.GetRMS();
206  }
207  Int_t inNeuron = 0;
208  Double_t diff = 0.;
209  if(fAnalysisTree) delete fAnalysisTree;
210  fAnalysisTree = new TTree("result","analysis");
212  fAnalysisTree->Branch("inNeuron",&inNeuron,"inNeuron/I");
213  fAnalysisTree->Branch("diff",&diff,"diff/D");
214  Int_t numOutNodes=GetNeurons(GetLayers());
215  Double_t *outVal=new Double_t[numOutNodes];
216  Double_t *trueVal=new Double_t[numOutNodes];
217 
218  delete fIOTree;
219  fIOTree=new TTree("MLP_iotree","MLP_iotree");
220  fIOTree->SetDirectory(0);
221  TString leaflist;
222  for (i=0; i<nn; i++)
223  leaflist+=Form("In%d/D:",i);
224  leaflist.Remove(leaflist.Length()-1);
225  fIOTree->Branch("In", params, leaflist);
226 
227  leaflist="";
228  for (i=0; i<numOutNodes; i++)
229  leaflist+=Form("Out%d/D:",i);
230  leaflist.Remove(leaflist.Length()-1);
231  fIOTree->Branch("Out", outVal, leaflist);
232 
233  leaflist="";
234  for (i=0; i<numOutNodes; i++)
235  leaflist+=Form("True%d/D:",i);
236  leaflist.Remove(leaflist.Length()-1);
237  fIOTree->Branch("True", trueVal, leaflist);
238  Double_t v1 = 0.;
239  Double_t v2 = 0.;
240  // Loop on the events in the test sample
241  for(j=0; j< nEvents; j++) {
242  fNetwork->GetEntry(test->GetEntry(j));
243  // Loop on the neurons to evaluate
244  for(k=0; k<GetNeurons(1); k++) {
245  params[k] = formulas[k]->EvalInstance(index[k]);
246  }
247  for(k=0; k<GetNeurons(GetLayers()); k++) {
248  outVal[k] = fNetwork->Evaluate(k,params);
249  trueVal[k] = ((TNeuron*)fNetwork->fLastLayer[k])->GetBranch();
250  }
251  fIOTree->Fill();
252 
253  // Loop on the input neurons
254  for (i = 0; i < GetNeurons(1); i++) {
255  inNeuron = i;
256  diff = 0;
257  // Loop on the neurons in the output layer
258  for(l=0; l<GetNeurons(GetLayers()); l++){
259  params[i] += shift*rms[i];
260  v1 = fNetwork->Evaluate(l,params);
261  params[i] -= 2*shift*rms[i];
262  v2 = fNetwork->Evaluate(l,params);
263  diff += (v1-v2)*(v1-v2);
264  // reset to original vealue
265  params[i] += shift*rms[i];
266  }
267  diff = TMath::Sqrt(diff);
268  fAnalysisTree->Fill();
269  }
270  }
271  delete[] params;
272  delete[] rms;
273  delete[] outVal;
274  delete[] trueVal;
275  delete[] index;
276  for(i=0; i<GetNeurons(1); i++) delete formulas[i]; delete [] formulas;
279 }
280 
281 ////////////////////////////////////////////////////////////////////////////////
282 /// Draws the distribution (on the test sample) of the
283 /// impact on the network output of a small variation of
284 /// the ith input.
285 
287 {
288  char sel[64];
289  snprintf(sel,64, "inNeuron==%d", i);
290  fAnalysisTree->Draw("diff", sel);
291 }
292 
293 ////////////////////////////////////////////////////////////////////////////////
294 /// Draws the distribution (on the test sample) of the
295 /// impact on the network output of a small variation of
296 /// each input.
297 /// DrawDInputs() draws something that approximates the distribution of the
298 /// derivative of the NN w.r.t. each input. That quantity is recognized as
299 /// one of the measures to determine key quantities in the network.
300 ///
301 /// What is done is to vary one input around its nominal value and to see
302 /// how the NN changes. This is done for each entry in the sample and produces
303 /// a distribution.
304 ///
305 /// What you can learn from that is:
306 /// - is variable a really useful, or is my network insensitive to it ?
307 /// - is there any risk of big systematic ? Is the network extremely sensitive
308 /// to small variations of any of my inputs ?
309 ///
310 /// As you might understand, this is to be considered with care and can serve
311 /// as input for an "educated guess" when optimizing the network.
312 
314 {
315  THStack* stack = new THStack("differences","differences (impact of variables on ANN)");
316  TLegend* legend = new TLegend(0.75,0.75,0.95,0.95);
317  TH1F* tmp = 0;
318  char var[64], sel[64];
319  for(Int_t i = 0; i < GetNeurons(1); i++) {
320  snprintf(var,64, "diff>>tmp%d", i);
321  snprintf(sel,64, "inNeuron==%d", i);
322  fAnalysisTree->Draw(var, sel, "goff");
323  tmp = (TH1F*)gDirectory->Get(Form("tmp%d",i));
324  tmp->SetDirectory(0);
325  tmp->SetLineColor(i+1);
326  stack->Add(tmp);
327  legend->AddEntry(tmp,GetInputNeuronTitle(i),"l");
328  }
329  stack->Draw("nostack");
330  legend->Draw();
331  gPad->SetLogy();
332 }
333 
334 ////////////////////////////////////////////////////////////////////////////////
335 /// Draws the distribution of the neural network (using ith neuron).
336 /// Two distributions are drawn, for events passing respectively the "signal"
337 /// and "background" cuts. Only the test sample is used.
338 
339 void TMLPAnalyzer::DrawNetwork(Int_t neuron, const char* signal, const char* bg)
340 {
341  TTree* data = fNetwork->fData;
342  TEventList* test = fNetwork->fTest;
343  TEventList* current = data->GetEventList();
344  data->SetEventList(test);
345  THStack* stack = new THStack("__NNout_TMLPA",Form("Neural net output (neuron %d)",neuron));
346  TH1F *bgh = new TH1F("__bgh_TMLPA", "NN output", 50, -0.5, 1.5);
347  TH1F *sigh = new TH1F("__sigh_TMLPA", "NN output", 50, -0.5, 1.5);
348  bgh->SetDirectory(0);
349  sigh->SetDirectory(0);
350  Int_t nEvents = 0;
351  Int_t j=0;
352  // build event lists for signal and background
353  TEventList* signal_list = new TEventList("__tmpSig_MLPA");
354  TEventList* bg_list = new TEventList("__tmpBkg_MLPA");
355  data->Draw(">>__tmpSig_MLPA",signal,"goff");
356  data->Draw(">>__tmpBkg_MLPA",bg,"goff");
357 
358  // fill the background
359  nEvents = bg_list->GetN();
360  for(j=0; j< nEvents; j++) {
361  bgh->Fill(fNetwork->Result(bg_list->GetEntry(j),neuron));
362  }
363  // fill the signal
364  nEvents = signal_list->GetN();
365  for(j=0; j< nEvents; j++) {
366  sigh->Fill(fNetwork->Result(signal_list->GetEntry(j),neuron));
367  }
368  // draws the result
369  bgh->SetLineColor(kBlue);
370  bgh->SetFillStyle(3008);
371  bgh->SetFillColor(kBlue);
372  sigh->SetLineColor(kRed);
373  sigh->SetFillStyle(3003);
374  sigh->SetFillColor(kRed);
375  bgh->SetStats(0);
376  sigh->SetStats(0);
377  stack->Add(bgh);
378  stack->Add(sigh);
379  TLegend *legend = new TLegend(.75, .80, .95, .95);
380  legend->AddEntry(bgh, "Background");
381  legend->AddEntry(sigh,"Signal");
382  stack->Draw("nostack");
383  legend->Draw();
384  // restore the default event list
385  data->SetEventList(current);
386  delete signal_list;
387  delete bg_list;
388 }
389 
390 ////////////////////////////////////////////////////////////////////////////////
391 /// Create a profile of the difference of the MLP output minus the
392 /// true value for a given output node outnode, vs the true value for
393 /// outnode, for all test data events. This method is mainly useful
394 /// when doing regression analysis with the MLP (i.e. not classification,
395 /// but continuous truth values).
396 /// The resulting TProfile histogram is returned.
397 /// It is not drawn if option "goff" is specified.
398 /// Options are passed to TProfile::Draw
399 
401  Option_t *option /*=""*/)
402 {
403  if (!fIOTree) GatherInformations();
404  TString pipehist=Form("MLP_truthdev_%d",outnode);
405  TString drawline;
406  drawline.Form("Out.Out%d-True.True%d:True.True%d>>",
407  outnode, outnode, outnode);
408  fIOTree->Draw(drawline+pipehist+"(20)", "", "goff prof");
409  TProfile* h=(TProfile*)gDirectory->Get(pipehist);
410  h->SetDirectory(0);
411  const char* title=GetOutputNeuronTitle(outnode);
412  if (title) {
413  h->SetTitle(Form("#Delta(output - truth) vs. truth for %s",
414  title));
415  h->GetXaxis()->SetTitle(title);
416  h->GetYaxis()->SetTitle(Form("#Delta(output - truth) for %s", title));
417  }
418  if (!strstr(option,"goff"))
419  h->Draw();
420  return h;
421 }
422 
423 ////////////////////////////////////////////////////////////////////////////////
424 /// Creates TProfiles of the difference of the MLP output minus the
425 /// true value vs the true value, one for each output, filled with the
426 /// test data events. This method is mainly useful when doing regression
427 /// analysis with the MLP (i.e. not classification, but continuous truth
428 /// values).
429 /// The returned THStack contains all the TProfiles. It is drawn unless
430 /// the option "goff" is specified.
431 /// Options are passed to TProfile::Draw.
432 
434 {
435  THStack *hs=new THStack("MLP_TruthDeviation",
436  "Deviation of MLP output from truth");
437 
438  // leg!=0 means we're drawing
439  TLegend *leg=0;
440  if (!option || !strstr(option,"goff"))
441  leg=new TLegend(.4,.85,.95,.95,"#Delta(output - truth) vs. truth for:");
442 
443  const char* xAxisTitle=0;
444 
445  // create profile for each input neuron,
446  // adding them into the THStack and the TLegend
447  for (Int_t outnode=0; outnode<GetNeurons(GetLayers()); outnode++) {
448  TProfile* h=DrawTruthDeviation(outnode, "goff");
449  h->SetLineColor(1+outnode);
450  hs->Add(h, option);
451  if (leg) leg->AddEntry(h,GetOutputNeuronTitle(outnode));
452  if (!outnode)
453  // Xaxis title is the same for all, extract it from the first one.
454  xAxisTitle=h->GetXaxis()->GetTitle();
455  }
456 
457  if (leg) {
458  hs->Draw("nostack");
459  leg->Draw();
460  // gotta draw before accessing the axes
461  hs->GetXaxis()->SetTitle(xAxisTitle);
462  hs->GetYaxis()->SetTitle("#Delta(output - truth)");
463  }
464 
465  return hs;
466 }
467 
468 ////////////////////////////////////////////////////////////////////////////////
469 /// Creates a profile of the difference of the MLP output outnode minus
470 /// the true value of outnode vs the input value innode, for all test
471 /// data events.
472 /// The resulting TProfile histogram is returned.
473 /// It is not drawn if option "goff" is specified.
474 /// Options are passed to TProfile::Draw
475 
477  Int_t outnode /*=0*/,
478  Option_t *option /*=""*/)
479 {
480  if (!fIOTree) GatherInformations();
481  TString pipehist=Form("MLP_truthdev_i%d_o%d", innode, outnode);
482  TString drawline;
483  drawline.Form("Out.Out%d-True.True%d:In.In%d>>",
484  outnode, outnode, innode);
485  fIOTree->Draw(drawline+pipehist+"(50)", "", "goff prof");
486  TProfile* h=(TProfile*)gROOT->FindObject(pipehist);
487  h->SetDirectory(0);
488  const char* titleInNeuron=GetInputNeuronTitle(innode);
489  const char* titleOutNeuron=GetOutputNeuronTitle(outnode);
490  h->SetTitle(Form("#Delta(output - truth) of %s vs. input %s",
491  titleOutNeuron, titleInNeuron));
492  h->GetXaxis()->SetTitle(Form("%s", titleInNeuron));
493  h->GetYaxis()->SetTitle(Form("#Delta(output - truth) for %s",
494  titleOutNeuron));
495  if (!strstr(option,"goff"))
496  h->Draw(option);
497  return h;
498 }
499 
500 ////////////////////////////////////////////////////////////////////////////////
501 /// Creates a profile of the difference of the MLP output outnode minus the
502 /// true value of outnode vs the input value, stacked for all inputs, for
503 /// all test data events.
504 /// The returned THStack contains all the TProfiles. It is drawn unless
505 /// the option "goff" is specified.
506 /// Options are passed to TProfile::Draw.
507 
509  Option_t *option /*=""*/)
510 {
511  TString sName;
512  sName.Form("MLP_TruthDeviationIO_%d", outnode);
513  const char* outputNodeTitle=GetOutputNeuronTitle(outnode);
514  THStack *hs=new THStack(sName,
515  Form("Deviation of MLP output %s from truth",
516  outputNodeTitle));
517 
518  // leg!=0 means we're drawing.
519  TLegend *leg=0;
520  if (!option || !strstr(option,"goff"))
521  leg=new TLegend(.4,.75,.95,.95,
522  Form("#Delta(output - truth) of %s vs. input for:",
523  outputNodeTitle));
524 
525  // create profile for each input neuron,
526  // adding them into the THStack and the TLegend
527  Int_t numInNodes=GetNeurons(1);
528  Int_t innode=0;
529  for (innode=0; innode<numInNodes; innode++) {
530  TProfile* h=DrawTruthDeviationInOut(innode, outnode, "goff");
531  h->SetLineColor(1+innode);
532  hs->Add(h, option);
533  if (leg) leg->AddEntry(h,h->GetXaxis()->GetTitle());
534  }
535 
536  if (leg) {
537  hs->Draw("nostack");
538  leg->Draw();
539  // gotta draw before accessing the axes
540  hs->GetXaxis()->SetTitle("Input value");
541  hs->GetYaxis()->SetTitle(Form("#Delta(output - truth) for %s",
542  outputNodeTitle));
543  }
544 
545  return hs;
546 }
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:851
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
const char * GetOutputNeuronTitle(Int_t out)
Returns the name of any neuron from the output layer.
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
const Double_t * v1
Definition: TArcBall.cxx:33
The Histogram stack class.
Definition: THStack.h:35
void CheckNetwork()
Gives some information about the network in the terminal.
TAxis * GetXaxis() const
Get x axis of the histogram used to draw the stack.
Definition: THStack.cxx:578
Ssiz_t Length() const
Definition: TString.h:390
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8266
const char Option_t
Definition: RtypesCore.h:62
const char * current
Definition: demos.C:12
TEventList * GetEventList() const
Definition: TTree.h:396
Definition: Rtypes.h:61
void DrawDInput(Int_t i)
Draws the distribution (on the test sample) of the impact on the network output of a small variation ...
virtual Int_t Fill()
Fill all branches.
Definition: TTree.cxx:4306
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
TH1 * h
Definition: legend2.C:5
Ssiz_t Index(const TString &str, Ssiz_t *len, Ssiz_t start=0) const
Find the first occurrence of the regexp in string and return the position, or -1 if there is no match...
Definition: TRegexp.cxx:208
Regular expression class.
Definition: TRegexp.h:35
#define gROOT
Definition: TROOT.h:344
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
int Int_t
Definition: RtypesCore.h:41
void GatherInformations()
Collect information about what is usefull in the network.
virtual void SetFillStyle(Style_t fstyle)
Definition: TAttFill.h:52
void DrawDInputs()
Draws the distribution (on the test sample) of the impact on the network output of a small variation ...
Profile Historam.
Definition: TProfile.h:34
Double_t Result(Int_t event, Int_t index=0) const
Computes the output for a given event.
Int_t GetLayers()
Returns the number of layers.
Int_t GetNeurons(Int_t layer)
Returns the number of neurons in given layer.
virtual void Draw(Option_t *chopt="")
Draw this multihist with its current attributes.
Definition: THStack.cxx:422
TLegend * legend
Definition: pirndm.C:35
void DrawNetwork(Int_t neuron, const char *signal, const char *bg)
Draws the distribution of the neural network (using ith neuron).
TProfile * DrawTruthDeviationInOut(Int_t innode, Int_t outnode=0, Option_t *option="")
Creates a profile of the difference of the MLP output outnode minus the true value of outnode vs the ...
virtual Int_t GetN() const
Definition: TEventList.h:58
const char * Data() const
Definition: TString.h:349
std::vector< std::vector< double > > Data
Int_t Atoi() const
Return integer value of string.
Definition: TString.cxx:1951
TMultiLayerPerceptron * fNetwork
Definition: TMLPAnalyzer.h:41
char * out
Definition: TBase64.cxx:29
Used to pass a selection expression to the Tree drawing routine.
Definition: TTreeFormula.h:64
virtual void SetEventList(TEventList *list)
This function transfroms the given TEventList into a TEntryList The new TEntryList is owned by the TT...
Definition: TTree.cxx:8164
TAxis * GetYaxis() const
Get x axis of the histogram used to draw the stack.
Definition: THStack.cxx:593
Double_t Evaluate(Int_t index, Double_t *params) const
Returns the Neural Net for a given set of input parameters #parameters must equal #input neurons...
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
const int nEvents
Definition: testRooFit.cxx:42
virtual Long64_t GetEntry(Int_t index) const
Return value of entry at index in the list.
Definition: TEventList.cxx:220
const char * GetInputNeuronTitle(Int_t in)
Returns the name of any neuron from the input layer.
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:7014
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2308
char * Form(const char *fmt,...)
A TEventList object is a list of selected events (entries) in a TTree.
Definition: TEventList.h:33
TLine * l
Definition: textangle.C:4
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
TAxis * GetYaxis()
Definition: TH1.h:320
const char * GetTitle() const
Returns title of object.
Definition: TAxis.h:133
Double_t GetRMS(Int_t axis=1) const
Definition: TH1.h:315
Int_t CountChar(Int_t c) const
Return number of times character c occurs in the string.
Definition: TString.cxx:430
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
TString & Remove(Ssiz_t pos)
Definition: TString.h:616
int Ssiz_t
Definition: RtypesCore.h:63
virtual void SetDirectory(TDirectory *dir)
Change the tree's directory.
Definition: TTree.cxx:8064
TTree * fAnalysisTree
Definition: TMLPAnalyzer.h:42
double Double_t
Definition: RtypesCore.h:55
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
leg
Definition: legend1.C:34
unsigned long ULong_t
Definition: RtypesCore.h:51
virtual void ResetBranchAddresses()
Tell all of our branches to drop their current objects and allocate new ones.
Definition: TTree.cxx:7266
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:360
virtual void Add(TH1 *h, Option_t *option="")
add a new histogram to the list Only 1-d and 2-d histograms currently supported.
Definition: THStack.cxx:336
T EvalInstance(Int_t i=0, const char *stringStack[]=0)
Evaluate this treeformula.
TEventList * fTest
EventList defining the events in the training dataset.
TString GetNeuronFormula(Int_t idx)
Returns the formula used as input for neuron (idx) in the first layer.
virtual Int_t Branch(TCollection *list, Int_t bufsize=32000, Int_t splitlevel=99, const char *name="")
Create one branch for each element in the collection.
Definition: TTree.cxx:1623
TProfile * DrawTruthDeviation(Int_t outnode=0, Option_t *option="")
Create a profile of the difference of the MLP output minus the true value for a given output node out...
#define gPad
Definition: TVirtualPad.h:288
A TTree object has a header with a name and a title.
Definition: TTree.h:98
#define gDirectory
Definition: TDirectory.h:221
THStack * DrawTruthDeviations(Option_t *option="")
Creates TProfiles of the difference of the MLP output minus the true value vs the true value...
Definition: Rtypes.h:61
void GetEntry(Int_t) const
Load an entry into the network.
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:6268
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
ClassImp(TMLPAnalyzer) TMLPAnalyzer
Destructor.
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
static void output(int code)
Definition: gifencode.c:226
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320
Ssiz_t First(char c) const
Find first occurrence of a character c.
Definition: TString.cxx:453
const char * cnt
Definition: TXMLSetup.cxx:75
TAxis * GetXaxis()
Definition: TH1.h:319
TTree * fIOTree
Definition: TMLPAnalyzer.h:43
THStack * DrawTruthDeviationInsOut(Int_t outnode=0, Option_t *option="")
Creates a profile of the difference of the MLP output outnode minus the true value of outnode vs the ...