ATLAS Offline Software
Loading...
Searching...
No Matches
TNetworkToHistoTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <TH1F.h>
6#include <TH2F.h>
7#include "TTrainedNetwork.h"
9#include <cmath>
10#include <vector>
11#include <iostream>
12
13using namespace std;
14
16
17std::vector<TH1*> TNetworkToHistoTool::fromTrainedNetworkToHisto(TTrainedNetwork* trainedNetwork) const
18{
19
20 std::vector<TH1*> outputHistos;
21
22 if (trainedNetwork->getActivationFunction()!=1)
23 {
24 cout << "ERROR: activation function is different from one. Only one is supported..." << endl;
25 return outputHistos;
26 }
27
28
29 Int_t nInput=trainedNetwork->getnInput();
30 vector<Int_t> nHiddenLayerSize=trainedNetwork->getnHiddenLayerSize();
31 Int_t nHidden=nHiddenLayerSize.size();
32
33 for (Int_t o=0;o<nHidden;++o)
34 {
35 cout << " Hidden lay: " << o << " size: " << nHiddenLayerSize[o];
36 }
37
38 Int_t nOutput=trainedNetwork->getnOutput();
39 cout << " Output size: " << nOutput << endl;
40
41 std::vector<TVectorD*> thresholdVectors=trainedNetwork->getThresholdVectors();
42 std::vector<TMatrixD*> weightMatrices=trainedNetwork->weightMatrices();
43
44 //LayersInfo
45
46 TH1F* histoLayersInfo=new TH1F("LayersInfo",
47 "LayersInfo",
48 nHidden+2,
49 0,
50 nHidden+2);
51
52 histoLayersInfo->SetBinContent(1,nInput);
53
54 for (Int_t i=0;i<nHidden;++i)
55 {
56 histoLayersInfo->SetBinContent(2+i,nHiddenLayerSize[i]);
57 }
58
59 histoLayersInfo->SetBinContent(2+nHidden,nOutput);
60
61 //underflow for linear output
62 if (trainedNetwork->getIfLinearOutput())
63 {
64 histoLayersInfo->SetBinContent(0,1);
65 }
66 //overflow for normalized output (Pott nodes)
67 if (trainedNetwork->getIfNormalizeOutput())
68 {
69 histoLayersInfo->SetBinContent(nHidden+3,1);
70 }
71
72 outputHistos.push_back(histoLayersInfo);
73
74
75 //ThresholdInfo
76 for (Int_t i=0;i<nHidden+1;++i)
77 {
78 TString threName("Layer");
79 threName+=i;
80 threName+="_thresholds";
81
82 Int_t layerSize=(i<nHidden)?nHiddenLayerSize[i]:nOutput;
83 Int_t previousLayerSize=(i==0)?nInput:nHiddenLayerSize[i-1];
84
85 TH1F* histoThreshLayer=new TH1F(threName,
86 threName,
87 layerSize,
88 0,
89 layerSize);
90
91 for (Int_t s=0;s<layerSize;s++)
92 {
93 histoThreshLayer->SetBinContent(s+1,thresholdVectors[i]->operator()(s));
94 }
95
96 TString weightsName("Layer");
97 weightsName+=i;
98 weightsName+="_weights";
99
100 outputHistos.push_back(histoThreshLayer);
101
102 TH2F* histoWeightsLayer=new TH2F(weightsName,
103 weightsName,
104 previousLayerSize,
105 0,
106 previousLayerSize,
107 layerSize,
108 0,
109 layerSize);
110
111 for (Int_t s=0;s<layerSize;s++)
112 {
113 for (Int_t p=0;p<previousLayerSize;++p)
114 {
115 histoWeightsLayer->SetBinContent(p+1,s+1,weightMatrices[i]->operator()(p,s));
116 }
117 }
118
119 outputHistos.push_back(histoWeightsLayer);
120
121 }
122
123
124 return outputHistos;
125
126}
127
128TH1* TNetworkToHistoTool::findHisto(TString nameOfHisto,
129 std::vector<TH1*> & inputHistos) const
130{
131
132 std::vector<TH1*>::const_iterator inputBegin=inputHistos.begin();
133 std::vector<TH1*>::const_iterator inputEnd=inputHistos.end();
134
135 for ( std::vector<TH1*>::const_iterator inputIter=inputBegin;inputIter!=inputEnd;++inputIter)
136 {
137 if ((*inputIter)->GetName()==nameOfHisto)
138 {
139 return (*inputIter);
140 }
141 }
142 return 0;
143}
144
145
146
147TTrainedNetwork* TNetworkToHistoTool::fromHistoToTrainedNetwork(std::vector<TH1*> & inputHistos) const
148{
149
150
151
152 TH1F* histoLayersInfo=dynamic_cast<TH1F*>(findHisto("LayersInfo",inputHistos));
153
154 if (histoLayersInfo==0)
155 {
156 cout << " Could not find LayersInfo histogram... Aborting " << endl;
157 return 0;
158 }
159
160 bool linearOutput=false;
161 bool normalizeOutput=false;
162
163
164 Int_t nHidden=histoLayersInfo->GetNbinsX()-2;
165 Int_t nInput=(Int_t)std::floor(histoLayersInfo->GetBinContent(1)+0.5);
166
167 vector<Int_t> nHiddenLayerSize;
168 for (Int_t i=0;i<nHidden;++i)
169 {
170 nHiddenLayerSize.push_back( (Int_t)std::floor(histoLayersInfo->GetBinContent(2+i)+0.5));
171 }
172
173 for (Int_t o=0;o<nHidden;++o)
174 {
175 cout << " Hidden lay: " << o << " size: " << nHiddenLayerSize[o];
176 }
177
178 Int_t nOutput=(Int_t)std::floor(histoLayersInfo->GetBinContent(2+nHidden)+0.5);
179 cout << " Output size: " << nOutput << endl;
180
181 if (histoLayersInfo->GetBinContent(0)>0.5)
182 {
183 linearOutput=true;
184 }
185 if (histoLayersInfo->GetBinContent(nHidden+3)>0.5)
186 {
187 normalizeOutput=true;
188 }
189
190 std::vector<TVectorD*> thresholdVectors;
191 std::vector<TMatrixD*> weightMatrices;
192
193
194 //Reconstruct thresholdInfo
195 for (Int_t i=0;i<nHidden+1;++i)
196 {
197 TString threName("Layer");
198 threName+=i;
199 threName+="_thresholds";
200
201 Int_t layerSize=(i<nHidden)?nHiddenLayerSize[i]:nOutput;
202 Int_t previousLayerSize=(i==0)?nInput:nHiddenLayerSize[i-1];
203
204 TVectorD* thresholdVector=new TVectorD(layerSize);
205 TMatrixD* weightMatrix=new TMatrixD(previousLayerSize,layerSize);
206
207 TH1F* histoThreshLayer=dynamic_cast<TH1F*>(findHisto(threName,inputHistos));
208 if (histoThreshLayer==0)
209 {
210 cout << " Could not find " << threName << " histogram... Aborting (mem leak also...)" << endl;
211 return 0;
212 }
213
214
215 for (Int_t s=0;s<layerSize;s++)
216 {
217 thresholdVector->operator()(s)=histoThreshLayer->GetBinContent(s+1);
218 }
219
220 TString weightsName("Layer");
221 weightsName+=i;
222 weightsName+="_weights";
223
224 TH2F* histoWeightsLayer=dynamic_cast<TH2F*>(findHisto(weightsName,inputHistos));
225 if (histoWeightsLayer==0)
226 {
227 cout << " Could not find " << weightsName << " histogram... Aborting (mem leak also...)" << endl;
228 return 0;
229 }
230
231 for (Int_t s=0;s<layerSize;s++)
232 {
233 for (Int_t p=0;p<previousLayerSize;++p)
234 {
235 weightMatrix->operator()(p,s)=histoWeightsLayer->GetBinContent(p+1,s+1);
236 }
237 }
238
239 thresholdVectors.push_back(thresholdVector);
240 weightMatrices.push_back(weightMatrix);
241
242 }
243
244
245 TTrainedNetwork* trainedNetwork=new TTrainedNetwork(nInput,
246 nHidden,
247 nOutput,
248 nHiddenLayerSize,
249 thresholdVectors,
250 weightMatrices,
251 1,
252 linearOutput,
253 normalizeOutput);
254 return trainedNetwork;
255
256}
257
258
259
ClassImp(xAOD::Experimental::RFileChecker) namespace xAOD
std::vector< TH1 * > fromTrainedNetworkToHisto(TTrainedNetwork *) const
TTrainedNetwork * fromHistoToTrainedNetwork(std::vector< TH1 * > &) const
TH1 * findHisto(TString nameOfHisto, std::vector< TH1 * > &inputHistos) const
TH2F(name, title, nxbins, bins_par2, bins_par3, bins_par4, bins_par5=None, bins_par6=None, path='', **kwargs)
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
STL namespace.