41 ClassImp(TMLPAnalyzer);
46 TMLPAnalyzer::~TMLPAnalyzer()
55 Int_t TMLPAnalyzer::GetLayers()
57 TString fStructure = fNetwork->GetStructure();
58 return fStructure.CountChar(
':')+1;
64 Int_t TMLPAnalyzer::GetNeurons(Int_t layer)
67 TString fStructure = fNetwork->GetStructure();
68 TString input = TString(fStructure(0, fStructure.First(
':')));
69 return input.CountChar(
',')+1;
71 else if(layer==GetLayers()) {
72 TString fStructure = fNetwork->GetStructure();
73 TString output = TString(fStructure(fStructure.Last(
':') + 1,
74 fStructure.Length() - fStructure.Last(
':')));
75 return output.CountChar(
',')+1;
79 TString fStructure = fNetwork->GetStructure();
80 TString hidden = TString(fStructure(fStructure.First(
':') + 1,
81 fStructure.Last(
':') - fStructure.First(
':') - 1));
83 Int_t end = hidden.Index(
":", beg + 1);
86 num = atoi(TString(hidden(beg, end - beg)).Data());
89 end = hidden.Index(
":", beg + 1);
90 if(layer==cnt)
return num;
92 num = atoi(TString(hidden(beg, hidden.Length() - beg)).Data());
94 if(layer==cnt)
return num;
103 TString TMLPAnalyzer::GetNeuronFormula(Int_t idx)
105 TString fStructure = fNetwork->GetStructure();
106 TString input = TString(fStructure(0, fStructure.First(
':')));
108 Int_t end = input.Index(
",", beg + 1);
112 brName = TString(input(beg, end - beg));
114 brName = brName(1,brName.Length()-1);
116 end = input.Index(
",", beg + 1);
117 if(cnt==idx)
return brName;
120 brName = TString(input(beg, input.Length() - beg));
122 brName = brName(1,brName.Length()-1);
129 const char* TMLPAnalyzer::GetInputNeuronTitle(Int_t in)
131 TNeuron* neuron=(TNeuron*)fNetwork->fFirstLayer[in];
132 return neuron ? neuron->GetName() :
"NO SUCH NEURON";
138 const char* TMLPAnalyzer::GetOutputNeuronTitle(Int_t out)
140 TNeuron* neuron=(TNeuron*)fNetwork->fLastLayer[out];
141 return neuron ? neuron->GetName() :
"NO SUCH NEURON";
147 void TMLPAnalyzer::CheckNetwork()
149 TString fStructure = fNetwork->GetStructure();
150 std::cout <<
"Network with structure: " << fStructure.Data() << std::endl;
151 std::cout <<
"inputs with low values in the differences plot may not be needed" << std::endl;
153 char var[64], sel[64];
154 for (Int_t i = 0; i < GetNeurons(1); i++) {
155 snprintf(var,64,
"diff>>tmp%d",i);
156 snprintf(sel,64,
"inNeuron==%d",i);
157 fAnalysisTree->Draw(var, sel,
"goff");
158 TH1F* tmp = (TH1F*)gDirectory->Get(Form(
"tmp%d",i));
160 std::cout << GetInputNeuronTitle(i)
161 <<
" -> " << tmp->GetMean()
162 <<
" +/- " << tmp->GetRMS() << std::endl;
171 void TMLPAnalyzer::GatherInformations()
173 Double_t shift = 0.1;
174 TTree* data = fNetwork->fData;
175 TEventList* test = fNetwork->fTest;
176 Int_t nEvents = test->GetN();
177 Int_t nn = GetNeurons(1);
178 Double_t* params =
new Double_t[nn];
179 Double_t* rms =
new Double_t[nn];
180 TTreeFormula** formulas =
new TTreeFormula*[nn];
181 Int_t* index =
new Int_t[nn];
183 TRegexp re(
"{[0-9]+}$");
184 Ssiz_t len = formula.Length();
186 Int_t i(0), j(0), k(0), l(0);
188 formula = GetNeuronFormula(i);
189 pos = re.Index(formula,&len);
190 if(pos==-1 || len<3) {
191 formulas[i] =
new TTreeFormula(Form(
"NF%lu",(ULong_t)
this),formula,data);
195 TString newformula(formula,pos);
196 TString val = formula(pos+1,len-2);
197 formulas[i] =
new TTreeFormula(Form(
"NF%lu",(ULong_t)
this),newformula,data);
198 formula = newformula;
199 index[i] = val.Atoi();
201 TH1D tmp(
"tmpb",
"tmpb", 1, -FLT_MAX, FLT_MAX);
202 data->Draw(Form(
"%s>>tmpb",formula.Data()),
"",
"goff");
203 rms[i] = tmp.GetRMS();
207 if(fAnalysisTree)
delete fAnalysisTree;
208 fAnalysisTree =
new TTree(
"result",
"analysis");
209 fAnalysisTree->SetDirectory(0);
210 fAnalysisTree->Branch(
"inNeuron",&inNeuron,
"inNeuron/I");
211 fAnalysisTree->Branch(
"diff",&diff,
"diff/D");
212 Int_t numOutNodes=GetNeurons(GetLayers());
213 Double_t *outVal=
new Double_t[numOutNodes];
214 Double_t *trueVal=
new Double_t[numOutNodes];
217 fIOTree=
new TTree(
"MLP_iotree",
"MLP_iotree");
218 fIOTree->SetDirectory(0);
221 leaflist+=Form(
"In%d/D:",i);
222 leaflist.Remove(leaflist.Length()-1);
223 fIOTree->Branch(
"In", params, leaflist);
226 for (i=0; i<numOutNodes; i++)
227 leaflist+=Form(
"Out%d/D:",i);
228 leaflist.Remove(leaflist.Length()-1);
229 fIOTree->Branch(
"Out", outVal, leaflist);
232 for (i=0; i<numOutNodes; i++)
233 leaflist+=Form(
"True%d/D:",i);
234 leaflist.Remove(leaflist.Length()-1);
235 fIOTree->Branch(
"True", trueVal, leaflist);
239 for(j=0; j< nEvents; j++) {
240 fNetwork->GetEntry(test->GetEntry(j));
242 for(k=0; k<GetNeurons(1); k++) {
243 params[k] = formulas[k]->EvalInstance(index[k]);
245 for(k=0; k<GetNeurons(GetLayers()); k++) {
246 outVal[k] = fNetwork->Evaluate(k,params);
247 trueVal[k] = ((TNeuron*)fNetwork->fLastLayer[k])->GetBranch();
252 for (i = 0; i < GetNeurons(1); i++) {
256 for(l=0; l<GetNeurons(GetLayers()); l++){
257 params[i] += shift*rms[i];
258 v1 = fNetwork->Evaluate(l,params);
259 params[i] -= 2*shift*rms[i];
260 v2 = fNetwork->Evaluate(l,params);
261 diff += (v1-v2)*(v1-v2);
263 params[i] += shift*rms[i];
265 diff = TMath::Sqrt(diff);
266 fAnalysisTree->Fill();
274 for(i=0; i<GetNeurons(1); i++)
delete formulas[i];
276 fAnalysisTree->ResetBranchAddresses();
277 fIOTree->ResetBranchAddresses();
285 void TMLPAnalyzer::DrawDInput(Int_t i)
288 snprintf(sel,64,
"inNeuron==%d", i);
289 fAnalysisTree->Draw(
"diff", sel);
312 void TMLPAnalyzer::DrawDInputs()
314 THStack* stack =
new THStack(
"differences",
"differences (impact of variables on ANN)");
315 TLegend* legend =
new TLegend(0.75,0.75,0.95,0.95);
317 char var[64], sel[64];
318 for(Int_t i = 0; i < GetNeurons(1); i++) {
319 snprintf(var,64,
"diff>>tmp%d", i);
320 snprintf(sel,64,
"inNeuron==%d", i);
321 fAnalysisTree->Draw(var, sel,
"goff");
322 tmp = (TH1F*)gDirectory->Get(Form(
"tmp%d",i));
323 tmp->SetDirectory(0);
324 tmp->SetLineColor(i+1);
326 legend->AddEntry(tmp,GetInputNeuronTitle(i),
"l");
328 stack->Draw(
"nostack");
338 void TMLPAnalyzer::DrawNetwork(Int_t neuron,
const char* signal,
const char* bg)
340 TTree* data = fNetwork->fData;
341 TEventList* test = fNetwork->fTest;
342 TEventList* current = data->GetEventList();
343 data->SetEventList(test);
344 THStack* stack =
new THStack(
"__NNout_TMLPA",Form(
"Neural net output (neuron %d)",neuron));
345 TH1F *bgh =
new TH1F(
"__bgh_TMLPA",
"NN output", 50, -0.5, 1.5);
346 TH1F *sigh =
new TH1F(
"__sigh_TMLPA",
"NN output", 50, -0.5, 1.5);
347 bgh->SetDirectory(0);
348 sigh->SetDirectory(0);
352 TEventList* signal_list =
new TEventList(
"__tmpSig_MLPA");
353 TEventList* bg_list =
new TEventList(
"__tmpBkg_MLPA");
354 data->Draw(
">>__tmpSig_MLPA",signal,
"goff");
355 data->Draw(
">>__tmpBkg_MLPA",bg,
"goff");
358 nEvents = bg_list->GetN();
359 for(j=0; j< nEvents; j++) {
360 bgh->Fill(fNetwork->Result(bg_list->GetEntry(j),neuron));
363 nEvents = signal_list->GetN();
364 for(j=0; j< nEvents; j++) {
365 sigh->Fill(fNetwork->Result(signal_list->GetEntry(j),neuron));
368 bgh->SetLineColor(kBlue);
369 bgh->SetFillStyle(3008);
370 bgh->SetFillColor(kBlue);
371 sigh->SetLineColor(kRed);
372 sigh->SetFillStyle(3003);
373 sigh->SetFillColor(kRed);
378 TLegend *legend =
new TLegend(.75, .80, .95, .95);
379 legend->AddEntry(bgh,
"Background");
380 legend->AddEntry(sigh,
"Signal");
381 stack->Draw(
"nostack");
384 data->SetEventList(current);
399 TProfile* TMLPAnalyzer::DrawTruthDeviation(Int_t outnode ,
402 if (!fIOTree) GatherInformations();
403 TString pipehist=Form(
"MLP_truthdev_%d",outnode);
405 drawline.Form(
"Out.Out%d-True.True%d:True.True%d>>",
406 outnode, outnode, outnode);
407 fIOTree->Draw(drawline+pipehist+
"(20)",
"",
"goff prof");
408 TProfile* h=(TProfile*)gDirectory->Get(pipehist);
410 const char* title=GetOutputNeuronTitle(outnode);
412 h->SetTitle(Form(
"#Delta(output - truth) vs. truth for %s",
414 h->GetXaxis()->SetTitle(title);
415 h->GetYaxis()->SetTitle(Form(
"#Delta(output - truth) for %s", title));
417 if (!strstr(option,
"goff"))
432 THStack* TMLPAnalyzer::DrawTruthDeviations(Option_t *option )
434 THStack *hs=
new THStack(
"MLP_TruthDeviation",
435 "Deviation of MLP output from truth");
439 if (!option || !strstr(option,
"goff"))
440 leg=
new TLegend(.4,.85,.95,.95,
"#Delta(output - truth) vs. truth for:");
442 const char* xAxisTitle=0;
446 for (Int_t outnode=0; outnode<GetNeurons(GetLayers()); outnode++) {
447 TProfile* h=DrawTruthDeviation(outnode,
"goff");
448 h->SetLineColor(1+outnode);
450 if (leg) leg->AddEntry(h,GetOutputNeuronTitle(outnode));
453 xAxisTitle=h->GetXaxis()->GetTitle();
460 hs->GetXaxis()->SetTitle(xAxisTitle);
461 hs->GetYaxis()->SetTitle(
"#Delta(output - truth)");
475 TProfile* TMLPAnalyzer::DrawTruthDeviationInOut(Int_t innode,
479 if (!fIOTree) GatherInformations();
480 TString pipehist=Form(
"MLP_truthdev_i%d_o%d", innode, outnode);
482 drawline.Form(
"Out.Out%d-True.True%d:In.In%d>>",
483 outnode, outnode, innode);
484 fIOTree->Draw(drawline+pipehist+
"(50)",
"",
"goff prof");
485 TProfile* h=(TProfile*)gROOT->FindObject(pipehist);
487 const char* titleInNeuron=GetInputNeuronTitle(innode);
488 const char* titleOutNeuron=GetOutputNeuronTitle(outnode);
489 h->SetTitle(Form(
"#Delta(output - truth) of %s vs. input %s",
490 titleOutNeuron, titleInNeuron));
491 h->GetXaxis()->SetTitle(Form(
"%s", titleInNeuron));
492 h->GetYaxis()->SetTitle(Form(
"#Delta(output - truth) for %s",
494 if (!strstr(option,
"goff"))
507 THStack* TMLPAnalyzer::DrawTruthDeviationInsOut(Int_t outnode ,
511 sName.Form(
"MLP_TruthDeviationIO_%d", outnode);
512 const char* outputNodeTitle=GetOutputNeuronTitle(outnode);
513 THStack *hs=
new THStack(sName,
514 Form(
"Deviation of MLP output %s from truth",
519 if (!option || !strstr(option,
"goff"))
520 leg=
new TLegend(.4,.75,.95,.95,
521 Form(
"#Delta(output - truth) of %s vs. input for:",
526 Int_t numInNodes=GetNeurons(1);
528 for (innode=0; innode<numInNodes; innode++) {
529 TProfile* h=DrawTruthDeviationInOut(innode, outnode,
"goff");
530 h->SetLineColor(1+innode);
532 if (leg) leg->AddEntry(h,h->GetXaxis()->GetTitle());
539 hs->GetXaxis()->SetTitle(
"Input value");
540 hs->GetYaxis()->SetTitle(Form(
"#Delta(output - truth) for %s",