ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
mlpHiggs.C
Go to the documentation of this file.
1 void mlpHiggs(Int_t ntrain=100) {
2 // Example of a Multi Layer Perceptron
3 // For a LEP search for invisible Higgs boson, a neural network
4 // was used to separate the signal from the background passing
5 // some selection cuts. Here is a simplified version of this network,
6 // taking into account only WW events.
7 //Author: Christophe Delaere
8 
9  if (!gROOT->GetClass("TMultiLayerPerceptron")) {
10  gSystem->Load("libMLP");
11  }
12 
13  // Prepare inputs
14  // The 2 trees are merged into one, and a "type" branch,
15  // equal to 1 for the signal and 0 for the background is added.
16  const char *fname = "mlpHiggs.root";
17  TFile *input = 0;
18  if (!gSystem->AccessPathName(fname)) {
19  input = TFile::Open(fname);
20  } else {
21  printf("accessing %s file from http://root.cern.ch/files\n",fname);
22  input = TFile::Open(Form("http://root.cern.ch/files/%s",fname));
23  }
24  if (!input) return;
25 
26  TTree *signal = (TTree *) input->Get("sig_filtered");
27  TTree *background = (TTree *) input->Get("bg_filtered");
28  TTree *simu = new TTree("MonteCarlo", "Filtered Monte Carlo Events");
29  Float_t ptsumf, qelep, nch, msumf, minvis, acopl, acolin;
30  Int_t type;
31  signal->SetBranchAddress("ptsumf", &ptsumf);
32  signal->SetBranchAddress("qelep", &qelep);
33  signal->SetBranchAddress("nch", &nch);
34  signal->SetBranchAddress("msumf", &msumf);
35  signal->SetBranchAddress("minvis", &minvis);
36  signal->SetBranchAddress("acopl", &acopl);
37  signal->SetBranchAddress("acolin", &acolin);
38  background->SetBranchAddress("ptsumf", &ptsumf);
39  background->SetBranchAddress("qelep", &qelep);
40  background->SetBranchAddress("nch", &nch);
41  background->SetBranchAddress("msumf", &msumf);
42  background->SetBranchAddress("minvis", &minvis);
43  background->SetBranchAddress("acopl", &acopl);
44  background->SetBranchAddress("acolin", &acolin);
45  simu->Branch("ptsumf", &ptsumf, "ptsumf/F");
46  simu->Branch("qelep", &qelep, "qelep/F");
47  simu->Branch("nch", &nch, "nch/F");
48  simu->Branch("msumf", &msumf, "msumf/F");
49  simu->Branch("minvis", &minvis, "minvis/F");
50  simu->Branch("acopl", &acopl, "acopl/F");
51  simu->Branch("acolin", &acolin, "acolin/F");
52  simu->Branch("type", &type, "type/I");
53  type = 1;
54  Int_t i;
55  for (i = 0; i < signal->GetEntries(); i++) {
56  signal->GetEntry(i);
57  simu->Fill();
58  }
59  type = 0;
60  for (i = 0; i < background->GetEntries(); i++) {
61  background->GetEntry(i);
62  simu->Fill();
63  }
64  // Build and train the NN ptsumf is used as a weight since we are primarly
65  // interested by high pt events.
66  // The datasets used here are the same as the default ones.
68  new TMultiLayerPerceptron("@msumf,@ptsumf,@acolin:5:3:type",
69  "ptsumf",simu,"Entry$%2","(Entry$+1)%2");
70  mlp->Train(ntrain, "text,graph,update=10");
71  mlp->Export("test","python");
72  // Use TMLPAnalyzer to see what it looks for
73  TCanvas* mlpa_canvas = new TCanvas("mlpa_canvas","Network analysis");
74  mlpa_canvas->Divide(2,2);
75  TMLPAnalyzer ana(mlp);
76  // Initialisation
77  ana.GatherInformations();
78  // output to the console
79  ana.CheckNetwork();
80  mlpa_canvas->cd(1);
81  // shows how each variable influences the network
82  ana.DrawDInputs();
83  mlpa_canvas->cd(2);
84  // shows the network structure
85  mlp->Draw();
86  mlpa_canvas->cd(3);
87  // draws the resulting network
88  ana.DrawNetwork(0,"type==1","type==0");
89  mlpa_canvas->cd(4);
90  // Use the NN to plot the results for each sample
91  // This will give approx. the same result as DrawNetwork.
92  // All entries are used, while DrawNetwork focuses on
93  // the test sample. Also the xaxis range is manually set.
94  TH1F *bg = new TH1F("bgh", "NN output", 50, -.5, 1.5);
95  TH1F *sig = new TH1F("sigh", "NN output", 50, -.5, 1.5);
96  bg->SetDirectory(0);
97  sig->SetDirectory(0);
98  Double_t params[3];
99  for (i = 0; i < background->GetEntries(); i++) {
100  background->GetEntry(i);
101  params[0] = msumf;
102  params[1] = ptsumf;
103  params[2] = acolin;
104  bg->Fill(mlp->Evaluate(0, params));
105  }
106  for (i = 0; i < signal->GetEntries(); i++) {
107  signal->GetEntry(i);
108  params[0] = msumf;
109  params[1] = ptsumf;
110  params[2] = acolin;
111  sig->Fill(mlp->Evaluate(0,params));
112  }
113  bg->SetLineColor(kBlue);
114  bg->SetFillStyle(3008); bg->SetFillColor(kBlue);
115  sig->SetLineColor(kRed);
116  sig->SetFillStyle(3003); sig->SetFillColor(kRed);
117  bg->SetStats(0);
118  sig->SetStats(0);
119  bg->Draw();
120  sig->Draw("same");
121  TLegend *legend = new TLegend(.75, .80, .95, .95);
122  legend->AddEntry(bg, "Background (WW)");
123  legend->AddEntry(sig, "Signal (Higgs)");
124  legend->Draw();
125  mlpa_canvas->cd(0);
126  delete input;
127 }
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1213
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
void CheckNetwork()
Gives some information about the network in the terminal.
float Float_t
Definition: RtypesCore.h:53
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
Definition: Rtypes.h:61
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
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
#define gROOT
Definition: TROOT.h:344
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all branches of entry and return total number of bytes read.
Definition: TTree.cxx:5144
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1766
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 ...
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).
void Export(Option_t *filename="NNfunction", Option_t *language="C++") const
Exports the NN as a function for any non-ROOT-dependant code Supported languages are: only C++ ...
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3851
void mlpHiggs(Int_t ntrain=100)
Definition: mlpHiggs.C:1
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=0)
Change branch address, dealing with clone trees properly.
Definition: TTree.cxx:7510
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
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
void Train(Int_t nEpoch, Option_t *option="text", Double_t minE=0)
Train the network.
virtual void Draw(Option_t *option="")
Draws the network structure.
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
char * Form(const char *fmt,...)
The Canvas class.
Definition: TCanvas.h:48
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
int type
Definition: TGX11.cxx:120
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
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
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1073
virtual Long64_t GetEntries() const
Definition: TTree.h:386
A TTree object has a header with a name and a title.
Definition: TTree.h:98
Definition: Rtypes.h:61
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320
TH1F * background
Definition: fithist.C:4