Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches

Detailed Description

View in nbviewer Open in SWAN Example of a Multi Layer Perceptron For a LEP search for invisible Higgs boson, a neural network was used to separate the signal from the background passing some selection cuts.

Here is a simplified version of this network, taking into account only WW events.

Training the Neural Network
Epoch: 0 learn=0.12648 test=0.125597
Epoch: 10 learn=0.0951981 test=0.0892184
Epoch: 20 learn=0.0930322 test=0.0884902
Epoch: 30 learn=0.0924559 test=0.0879593
Epoch: 40 learn=0.0912523 test=0.0870583
Epoch: 50 learn=0.0909583 test=0.0870296
Epoch: 60 learn=0.0906589 test=0.0864934
Epoch: 70 learn=0.0904439 test=0.0865701
Epoch: 80 learn=0.0903373 test=0.0863813
Epoch: 90 learn=0.0896677 test=0.0853879
Epoch: 99 learn=0.0889658 test=0.0850847
Training done.
test.py created.
Network with structure: @msumf,@ptsumf,@acolin:5:3:type
inputs with low values in the differences plot may not be needed
@msumf -> 0.0150092 +/- 0.0192632
@ptsumf -> 0.0285032 +/- 0.041759
@acolin -> 0.0262982 +/- 0.0496902
void mlpHiggs(Int_t ntrain=100) {
const char *fname = "mlpHiggs.root";
TFile *input = 0;
if (!gSystem->AccessPathName(fname)) {
input = TFile::Open(fname);
} else if (!gSystem->AccessPathName(Form("%s/legacy/mlp/%s", TROOT::GetTutorialDir().Data(), fname))) {
input = TFile::Open(Form("%s/legacy/mlp/%s", TROOT::GetTutorialDir().Data(), fname));
} else {
printf("accessing %s file from http://root.cern.ch/files\n",fname);
input = TFile::Open(Form("http://root.cern.ch/files/%s",fname));
}
if (!input) return;
TTree *sig_filtered = (TTree *) input->Get("sig_filtered");
TTree *bg_filtered = (TTree *) input->Get("bg_filtered");
TTree *simu = new TTree("MonteCarlo", "Filtered Monte Carlo Events");
Float_t ptsumf, qelep, nch, msumf, minvis, acopl, acolin;
sig_filtered->SetBranchAddress("ptsumf", &ptsumf);
sig_filtered->SetBranchAddress("qelep", &qelep);
sig_filtered->SetBranchAddress("nch", &nch);
sig_filtered->SetBranchAddress("msumf", &msumf);
sig_filtered->SetBranchAddress("minvis", &minvis);
sig_filtered->SetBranchAddress("acopl", &acopl);
sig_filtered->SetBranchAddress("acolin", &acolin);
bg_filtered->SetBranchAddress("ptsumf", &ptsumf);
bg_filtered->SetBranchAddress("qelep", &qelep);
bg_filtered->SetBranchAddress("nch", &nch);
bg_filtered->SetBranchAddress("msumf", &msumf);
bg_filtered->SetBranchAddress("minvis", &minvis);
bg_filtered->SetBranchAddress("acopl", &acopl);
bg_filtered->SetBranchAddress("acolin", &acolin);
simu->Branch("ptsumf", &ptsumf, "ptsumf/F");
simu->Branch("qelep", &qelep, "qelep/F");
simu->Branch("nch", &nch, "nch/F");
simu->Branch("msumf", &msumf, "msumf/F");
simu->Branch("minvis", &minvis, "minvis/F");
simu->Branch("acopl", &acopl, "acopl/F");
simu->Branch("acolin", &acolin, "acolin/F");
simu->Branch("type", &type, "type/I");
type = 1;
Int_t i;
for (i = 0; i < sig_filtered->GetEntries(); i++) {
sig_filtered->GetEntry(i);
simu->Fill();
}
type = 0;
for (i = 0; i < bg_filtered->GetEntries(); i++) {
bg_filtered->GetEntry(i);
simu->Fill();
}
// Build and train the NN ptsumf is used as a weight since we are primarily
// interested by high pt events.
// The datasets used here are the same as the default ones.
new TMultiLayerPerceptron("@msumf,@ptsumf,@acolin:5:3:type",
"ptsumf",simu,"Entry$%2","(Entry$+1)%2");
mlp->Train(ntrain, "text,graph,update=10");
mlp->Export("test","python");
// Use TMLPAnalyzer to see what it looks for
TCanvas* mlpa_canvas = new TCanvas("mlpa_canvas","Network analysis");
mlpa_canvas->Divide(2,2);
TMLPAnalyzer ana(mlp);
// Initialisation
ana.GatherInformations();
// output to the console
ana.CheckNetwork();
mlpa_canvas->cd(1);
// shows how each variable influences the network
ana.DrawDInputs();
mlpa_canvas->cd(2);
// shows the network structure
mlp->Draw();
mlpa_canvas->cd(3);
// draws the resulting network
ana.DrawNetwork(0,"type==1","type==0");
mlpa_canvas->cd(4);
// Use the NN to plot the results for each sample
// This will give approx. the same result as DrawNetwork.
// All entries are used, while DrawNetwork focuses on
// the test sample. Also the xaxis range is manually set.
TH1F *bg = new TH1F("bgh", "NN output", 50, -.5, 1.5);
TH1F *sig = new TH1F("sigh", "NN output", 50, -.5, 1.5);
bg->SetDirectory(0);
sig->SetDirectory(0);
Double_t params[3];
for (i = 0; i < bg_filtered->GetEntries(); i++) {
bg_filtered->GetEntry(i);
params[0] = msumf;
params[1] = ptsumf;
params[2] = acolin;
bg->Fill(mlp->Evaluate(0, params));
}
for (i = 0; i < sig_filtered->GetEntries(); i++) {
sig_filtered->GetEntry(i);
params[0] = msumf;
params[1] = ptsumf;
params[2] = acolin;
sig->Fill(mlp->Evaluate(0,params));
}
bg->SetFillStyle(3008); bg->SetFillColor(kBlue);
sig->SetFillStyle(3003); sig->SetFillColor(kRed);
bg->SetStats(0);
sig->SetStats(0);
bg->Draw();
sig->Draw("same");
TLegend *legend = new TLegend(.75, .80, .95, .95);
legend->AddEntry(bg, "Background (WW)");
legend->AddEntry(sig, "Signal (Higgs)");
legend->Draw();
mlpa_canvas->cd(0);
delete input;
}
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
int type
Definition TGX11.cxx:121
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
Definition TSystem.h:559
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:708
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:54
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:3997
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
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:8777
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3350
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3073
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8830
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:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition TLegend.cxx:423
This utility class contains a set of tests usefull 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.
virtual void Draw(Option_t *option="")
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:1177
static const TString & GetTutorialDir()
Get the tutorials directory in the installation. Static utility function.
Definition TROOT.cxx:3032
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:1294
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Int_t Fill()
Fill all branches.
Definition TTree.cxx:4590
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=0)
Change branch address, dealing with clone trees properly.
Definition TTree.cxx:8349
virtual Long64_t GetEntries() const
Definition TTree.h:460
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:350
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:5618
Author
Christophe Delaere

Definition in file mlpHiggs.C.