Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
mlpHiggs.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_mlp
3/// \notebook
4/// Example of a Multi Layer Perceptron
5/// For a LEP search for invisible Higgs boson, a neural network
6/// was used to separate the signal from the background passing
7/// some selection cuts. Here is a simplified version of this network,
8/// taking into account only WW events.
9///
10/// \macro_image
11/// \macro_output
12/// \macro_code
13///
14/// \author Christophe Delaere
15
16void mlpHiggs(Int_t ntrain=100) {
17 const char *fname = "mlpHiggs.root";
18 TFile *input = 0;
19 if (!gSystem->AccessPathName(fname)) {
20 input = TFile::Open(fname);
21 } else if (!gSystem->AccessPathName(Form("%s/legacy/mlp/%s", TROOT::GetTutorialDir().Data(), fname))) {
22 input = TFile::Open(Form("%s/legacy/mlp/%s", TROOT::GetTutorialDir().Data(), fname));
23 } else {
24 printf("accessing %s file from http://root.cern/files\n",fname);
25 input = TFile::Open(Form("http://root.cern/files/%s",fname));
26 }
27 if (!input) return;
28
29 TTree *sig_filtered = (TTree *) input->Get("sig_filtered");
30 TTree *bg_filtered = (TTree *) input->Get("bg_filtered");
31 TTree *simu = new TTree("MonteCarlo", "Filtered Monte Carlo Events");
32 Float_t ptsumf, qelep, nch, msumf, minvis, acopl, acolin;
33 Int_t type;
34 sig_filtered->SetBranchAddress("ptsumf", &ptsumf);
35 sig_filtered->SetBranchAddress("qelep", &qelep);
36 sig_filtered->SetBranchAddress("nch", &nch);
37 sig_filtered->SetBranchAddress("msumf", &msumf);
38 sig_filtered->SetBranchAddress("minvis", &minvis);
39 sig_filtered->SetBranchAddress("acopl", &acopl);
40 sig_filtered->SetBranchAddress("acolin", &acolin);
41 bg_filtered->SetBranchAddress("ptsumf", &ptsumf);
42 bg_filtered->SetBranchAddress("qelep", &qelep);
43 bg_filtered->SetBranchAddress("nch", &nch);
44 bg_filtered->SetBranchAddress("msumf", &msumf);
45 bg_filtered->SetBranchAddress("minvis", &minvis);
46 bg_filtered->SetBranchAddress("acopl", &acopl);
47 bg_filtered->SetBranchAddress("acolin", &acolin);
48 simu->Branch("ptsumf", &ptsumf, "ptsumf/F");
49 simu->Branch("qelep", &qelep, "qelep/F");
50 simu->Branch("nch", &nch, "nch/F");
51 simu->Branch("msumf", &msumf, "msumf/F");
52 simu->Branch("minvis", &minvis, "minvis/F");
53 simu->Branch("acopl", &acopl, "acopl/F");
54 simu->Branch("acolin", &acolin, "acolin/F");
55 simu->Branch("type", &type, "type/I");
56 type = 1;
57 Int_t i;
58 for (i = 0; i < sig_filtered->GetEntries(); i++) {
59 sig_filtered->GetEntry(i);
60 simu->Fill();
61 }
62 type = 0;
63 for (i = 0; i < bg_filtered->GetEntries(); i++) {
64 bg_filtered->GetEntry(i);
65 simu->Fill();
66 }
67 // Build and train the NN ptsumf is used as a weight since we are primarily
68 // interested by high pt events.
69 // The datasets used here are the same as the default ones.
71 new TMultiLayerPerceptron("@msumf,@ptsumf,@acolin:5:3:type",
72 "ptsumf",simu,"Entry$%2","(Entry$+1)%2");
73 mlp->Train(ntrain, "text,graph,update=10");
74 mlp->Export("test","python");
75 // Use TMLPAnalyzer to see what it looks for
76 TCanvas* mlpa_canvas = new TCanvas("mlpa_canvas","Network analysis");
77 mlpa_canvas->Divide(2,2);
78 TMLPAnalyzer ana(mlp);
79 // Initialisation
80 ana.GatherInformations();
81 // output to the console
82 ana.CheckNetwork();
83 mlpa_canvas->cd(1);
84 // shows how each variable influences the network
85 ana.DrawDInputs();
86 mlpa_canvas->cd(2);
87 // shows the network structure
88 mlp->Draw();
89 mlpa_canvas->cd(3);
90 // draws the resulting network
91 ana.DrawNetwork(0,"type==1","type==0");
92 mlpa_canvas->cd(4);
93 // Use the NN to plot the results for each sample
94 // This will give approx. the same result as DrawNetwork.
95 // All entries are used, while DrawNetwork focuses on
96 // the test sample. Also the xaxis range is manually set.
97 TH1F *bg = new TH1F("bgh", "NN output", 50, -.5, 1.5);
98 TH1F *sig = new TH1F("sigh", "NN output", 50, -.5, 1.5);
99 bg->SetDirectory(0);
100 sig->SetDirectory(0);
101 Double_t params[3];
102 for (i = 0; i < bg_filtered->GetEntries(); i++) {
103 bg_filtered->GetEntry(i);
104 params[0] = msumf;
105 params[1] = ptsumf;
106 params[2] = acolin;
107 bg->Fill(mlp->Evaluate(0, params));
108 }
109 for (i = 0; i < sig_filtered->GetEntries(); i++) {
110 sig_filtered->GetEntry(i);
111 params[0] = msumf;
112 params[1] = ptsumf;
113 params[2] = acolin;
114 sig->Fill(mlp->Evaluate(0,params));
115 }
116 bg->SetLineColor(kBlue);
117 bg->SetFillStyle(3008); bg->SetFillColor(kBlue);
118 sig->SetLineColor(kRed);
119 sig->SetFillStyle(3003); sig->SetFillColor(kRed);
120 bg->SetStats(0);
121 sig->SetStats(0);
122 bg->Draw();
123 sig->Draw("same");
124 TLegend *legend = new TLegend(.75, .80, .95, .95);
125 legend->AddEntry(bg, "Background (WW)");
126 legend->AddEntry(sig, "Signal (Higgs)");
127 legend->Draw();
128 mlpa_canvas->cd(0);
129 delete input;
130}
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void input
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
R__EXTERN TSystem * gSystem
Definition TSystem.h:561
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:39
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:719
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4089
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:622
virtual void SetDirectory(TDirectory *dir)
By default, when a histogram is created, it is added to the list of histogram objects in the current ...
Definition TH1.cxx:8937
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3344
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8990
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition TLegend.cxx:320
void Draw(Option_t *option="") override
Draw this legend with its current attributes.
Definition TLegend.cxx:425
This utility class contains a set of tests useful when developing a neural network.
This class describes a neural network.
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.
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++ ,...
void Train(Int_t nEpoch, Option_t *option="text", Double_t minE=0)
Train the network.
void Draw(Option_t *option="") override
Draws the network structure.
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
Definition TPad.cxx:1249
static const TString & GetTutorialDir()
Get the tutorials directory in the installation. Static utility function.
Definition TROOT.cxx:3119
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:1296
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Int_t Fill()
Fill all branches.
Definition TTree.cxx:4603
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
Read all branches of entry and return total number of bytes read.
Definition TTree.cxx:5638
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=nullptr)
Change branch address, dealing with clone trees properly.
Definition TTree.cxx:8385
virtual Long64_t GetEntries() const
Definition TTree.h:463
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition TTree.h:353