00001 void mlpHiggs(Int_t ntrain=100) {
00002
00003
00004
00005
00006
00007
00008
00009 if (!gROOT->GetClass("TMultiLayerPerceptron")) {
00010 gSystem->Load("libMLP");
00011 }
00012
00013
00014
00015
00016 const char *fname = "mlpHiggs.root";
00017 TFile *input = 0;
00018 if (!gSystem->AccessPathName(fname)) {
00019 input = TFile::Open(fname);
00020 } else {
00021 printf("accessing %s file from http://root.cern.ch/files\n",fname);
00022 input = TFile::Open(Form("http://root.cern.ch/files/%s",fname));
00023 }
00024 if (!input) return;
00025
00026 TTree *signal = (TTree *) input->Get("sig_filtered");
00027 TTree *background = (TTree *) input->Get("bg_filtered");
00028 TTree *simu = new TTree("MonteCarlo", "Filtered Monte Carlo Events");
00029 Float_t ptsumf, qelep, nch, msumf, minvis, acopl, acolin;
00030 Int_t type;
00031 signal->SetBranchAddress("ptsumf", &ptsumf);
00032 signal->SetBranchAddress("qelep", &qelep);
00033 signal->SetBranchAddress("nch", &nch);
00034 signal->SetBranchAddress("msumf", &msumf);
00035 signal->SetBranchAddress("minvis", &minvis);
00036 signal->SetBranchAddress("acopl", &acopl);
00037 signal->SetBranchAddress("acolin", &acolin);
00038 background->SetBranchAddress("ptsumf", &ptsumf);
00039 background->SetBranchAddress("qelep", &qelep);
00040 background->SetBranchAddress("nch", &nch);
00041 background->SetBranchAddress("msumf", &msumf);
00042 background->SetBranchAddress("minvis", &minvis);
00043 background->SetBranchAddress("acopl", &acopl);
00044 background->SetBranchAddress("acolin", &acolin);
00045 simu->Branch("ptsumf", &ptsumf, "ptsumf/F");
00046 simu->Branch("qelep", &qelep, "qelep/F");
00047 simu->Branch("nch", &nch, "nch/F");
00048 simu->Branch("msumf", &msumf, "msumf/F");
00049 simu->Branch("minvis", &minvis, "minvis/F");
00050 simu->Branch("acopl", &acopl, "acopl/F");
00051 simu->Branch("acolin", &acolin, "acolin/F");
00052 simu->Branch("type", &type, "type/I");
00053 type = 1;
00054 Int_t i;
00055 for (i = 0; i < signal->GetEntries(); i++) {
00056 signal->GetEntry(i);
00057 simu->Fill();
00058 }
00059 type = 0;
00060 for (i = 0; i < background->GetEntries(); i++) {
00061 background->GetEntry(i);
00062 simu->Fill();
00063 }
00064
00065
00066
00067 TMultiLayerPerceptron *mlp =
00068 new TMultiLayerPerceptron("@msumf,@ptsumf,@acolin:5:3:type",
00069 "ptsumf",simu,"Entry$%2","(Entry$+1)%2");
00070 mlp->Train(ntrain, "text,graph,update=10");
00071 mlp->Export("test","python");
00072
00073 TCanvas* mlpa_canvas = new TCanvas("mlpa_canvas","Network analysis");
00074 mlpa_canvas->Divide(2,2);
00075 TMLPAnalyzer ana(mlp);
00076
00077 ana.GatherInformations();
00078
00079 ana.CheckNetwork();
00080 mlpa_canvas->cd(1);
00081
00082 ana.DrawDInputs();
00083 mlpa_canvas->cd(2);
00084
00085 mlp->Draw();
00086 mlpa_canvas->cd(3);
00087
00088 ana.DrawNetwork(0,"type==1","type==0");
00089 mlpa_canvas->cd(4);
00090
00091
00092
00093
00094 TH1F *bg = new TH1F("bgh", "NN output", 50, -.5, 1.5);
00095 TH1F *sig = new TH1F("sigh", "NN output", 50, -.5, 1.5);
00096 bg->SetDirectory(0);
00097 sig->SetDirectory(0);
00098 Double_t params[4];
00099 for (i = 0; i < background->GetEntries(); i++) {
00100 background->GetEntry(i);
00101 params[0] = msumf;
00102 params[1] = ptsumf;
00103 params[2] = acolin;
00104 params[3] = acopl;
00105 bg->Fill(mlp->Evaluate(0, params));
00106 }
00107 for (i = 0; i < signal->GetEntries(); i++) {
00108 signal->GetEntry(i);
00109 params[0] = msumf;
00110 params[1] = ptsumf;
00111 params[2] = acolin;
00112 params[3] = acopl;
00113 sig->Fill(mlp->Evaluate(0,params));
00114 }
00115 bg->SetLineColor(kBlue);
00116 bg->SetFillStyle(3008); bg->SetFillColor(kBlue);
00117 sig->SetLineColor(kRed);
00118 sig->SetFillStyle(3003); sig->SetFillColor(kRed);
00119 bg->SetStats(0);
00120 sig->SetStats(0);
00121 bg->Draw();
00122 sig->Draw("same");
00123 TLegend *legend = new TLegend(.75, .80, .95, .95);
00124 legend->AddEntry(bg, "Background (WW)");
00125 legend->AddEntry(sig, "Signal (Higgs)");
00126 legend->Draw();
00127 mlpa_canvas->cd(0);
00128 delete input;
00129 }