Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMultiLayerPerceptron.cxx
Go to the documentation of this file.
1// @(#)root/mlp:$Id$
2// Author: Christophe.Delaere@cern.ch 20/07/03
3
4/*************************************************************************
5 * Copyright (C) 1995-2003, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TMultiLayerPerceptron
13
14
15This class describes a neural network.
16There are facilities to train the network and use the output.
17
18The input layer is made of inactive neurons (returning the
19optionally normalized input) and output neurons are linear.
20The type of hidden neurons is free, the default being sigmoids.
21(One should still try to pass normalized inputs, e.g. between [0.,1])
22
23The basic input is a TTree and two (training and test) TEventLists.
24Input and output neurons are assigned a value computed for each event
25with the same possibilities as for TTree::Draw().
26Events may be weighted individually or via TTree::SetWeight().
276 learning methods are available: kStochastic, kBatch,
28kSteepestDescent, kRibierePolak, kFletcherReeves and kBFGS.
29
30This implementation, written by C. Delaere, is *inspired* from
31the mlpfit package from J.Schwindling et al. with some extensions:
32
33 - the algorithms are globally the same
34 - in TMultilayerPerceptron, there is no limitation on the number of
35 layers/neurons, while MLPFIT was limited to 2 hidden layers
36 - TMultilayerPerceptron allows you to save the network in a root file, and
37 provides more export functionalities
38 - TMultilayerPerceptron gives more flexibility regarding the normalization of
39 inputs/outputs
40 - TMultilayerPerceptron provides, thanks to Andrea Bocci, the possibility to
41 use cross-entropy errors, which allows to train a network for pattern
42 classification based on Bayesian posterior probability.
43
44### Introduction
45
46Neural Networks are more and more used in various fields for data
47analysis and classification, both for research and commercial
48institutions. Some randomly chosen examples are:
49
50 - image analysis
51 - financial movements predictions and analysis
52 - sales forecast and product shipping optimisation
53 - in particles physics: mainly for classification tasks (signal
54 over background discrimination)
55
56More than 50% of neural networks are multilayer perceptrons. This
57implementation of multilayer perceptrons is inspired from the
58<A HREF="http://schwind.home.cern.ch/schwind/MLPfit.html">MLPfit
59package</A> originally written by Jerome Schwindling. MLPfit remains
60one of the fastest tool for neural networks studies, and this ROOT
61add-on will not try to compete on that. A clear and flexible Object
62Oriented implementation has been chosen over a faster but more
63difficult to maintain code. Nevertheless, the time penalty does not
64exceed a factor 2.
65
66### The MLP
67
68The multilayer perceptron is a simple feed-forward network with
69the following structure:
70
71\image html mlp.png
72
73It is made of neurons characterized by a bias and weighted links
74between them (let's call those links synapses). The input neurons
75receive the inputs, normalize them and forward them to the first
76hidden layer.
77
78Each neuron in any subsequent layer first computes a linear
79combination of the outputs of the previous layer. The output of the
80neuron is then function of that combination with <I>f</I> being
81linear for output neurons or a sigmoid for hidden layers. This is
82useful because of two theorems:
83
84 1. A linear combination of sigmoids can approximate any
85 continuous function.
86 2. Trained with output = 1 for the signal and 0 for the
87 background, the approximated function of inputs X is the probability
88 of signal, knowing X.
89
90### Learning methods
91
92The aim of all learning methods is to minimize the total error on
93a set of weighted examples. The error is defined as the sum in
94quadrature, divided by two, of the error on each individual output
95neuron.
96In all methods implemented, one needs to compute
97the first derivative of that error with respect to the weights.
98Exploiting the well-known properties of the derivative, especially the
99derivative of compound functions, one can write:
100
101 - for a neuron: product of the local derivative with the
102 weighted sum on the outputs of the derivatives.
103 - for a synapse: product of the input with the local derivative
104 of the output neuron.
105
106This computation is called back-propagation of the errors. A
107loop over all examples is called an epoch.
108Six learning methods are implemented.
109
110#### Stochastic minimization:
111
112is the most trivial learning method. This is the Robbins-Monro
113stochastic approximation applied to multilayer perceptrons. The
114weights are updated after each example according to the formula:
115\f$w_{ij}(t+1) = w_{ij}(t) + \Delta w_{ij}(t)\f$
116
117with
118
119\f$\Delta w_{ij}(t) = - \eta(d e_p / d w_{ij} + \delta) + \epsilon \Delta w_{ij}(t-1)\f$
120
121The parameters for this method are Eta, EtaDecay, Delta and
122Epsilon.
123
124#### Steepest descent with fixed step size (batch learning):
125
126It is the same as the stochastic
127minimization, but the weights are updated after considering all the
128examples, with the total derivative dEdw. The parameters for this
129method are Eta, EtaDecay, Delta and Epsilon.
130
131#### Steepest descent algorithm:
132
133Weights are set to the minimum along the line defined by the gradient. The
134only parameter for this method is Tau. Lower tau = higher precision =
135slower search. A value Tau = 3 seems reasonable.
136
137#### Conjugate gradients with the Polak-Ribiere updating formula:
138
139Weights are set to the minimum along the line defined by the conjugate gradient.
140Parameters are Tau and Reset, which defines the epochs where the direction is
141reset to the steepest descent.
142
143#### Conjugate gradients with the Fletcher-Reeves updating formula:
144
145Weights are set to the minimum along the line defined by the conjugate gradient. Parameters
146are Tau and Reset, which defines the epochs where the direction is
147reset to the steepest descent.
148
149#### Broyden, Fletcher, Goldfarb, Shanno (BFGS) method:
150
151 Implies the computation of a NxN matrix
152computation, but seems more powerful at least for less than 300
153weights. Parameters are Tau and Reset, which defines the epochs where
154the direction is reset to the steepest descent.
155
156### How to use it...
157
158TMLP is build from 3 classes: TNeuron, TSynapse and
159TMultiLayerPerceptron. Only TMultiLayerPerceptron should be used
160explicitly by the user.
161
162TMultiLayerPerceptron will take examples from a TTree
163given in the constructor. The network is described by a simple
164string: The input/output layers are defined by giving the expression for
165each neuron, separated by comas. Hidden layers are just described
166by the number of neurons. The layers are separated by colons.
167In addition, input/output layer formulas can be preceded by '@' (e.g "@out")
168if one wants to also normalize the data from the TTree.
169Input and outputs are taken from the TTree given as second argument.
170Expressions are evaluated as for TTree::Draw(), arrays are expended in
171distinct neurons, one for each index.
172This can only be done for fixed-size arrays.
173If the formula ends with "!", softmax functions are used for the output layer.
174One defines the training and test datasets by TEventLists.
175
176Example:
177~~~ {.cpp}
178TMultiLayerPerceptron("x,y:10:5:f",inputTree);
179~~~
180
181Both the TTree and the TEventLists can be defined in
182the constructor, or later with the suited setter method. The lists
183used for training and test can be defined either explicitly, or via
184a string containing the formula to be used to define them, exactly as
185for a TCut.
186
187The learning method is defined using the TMultiLayerPerceptron::SetLearningMethod() .
188Learning methods are :
189
190 - TMultiLayerPerceptron::kStochastic,
191 - TMultiLayerPerceptron::kBatch,
192 - TMultiLayerPerceptron::kSteepestDescent,
193 - TMultiLayerPerceptron::kRibierePolak,
194 - TMultiLayerPerceptron::kFletcherReeves,
195 - TMultiLayerPerceptron::kBFGS
196
197A weight can be assigned to events, either in the constructor, either
198with TMultiLayerPerceptron::SetEventWeight(). In addition, the TTree weight
199is taken into account.
200
201Finally, one starts the training with
202TMultiLayerPerceptron::Train(Int_t nepoch, Option_t* options). The
203first argument is the number of epochs while option is a string that
204can contain: "text" (simple text output) , "graph"
205(evoluting graphical training curves), "update=X" (step for
206the text/graph output update) or "+" (will skip the
207randomisation and start from the previous values). All combinations
208are available.
209
210Example:
211~~~ {.cpp}
212net.Train(100,"text, graph, update=10");
213~~~
214
215When the neural net is trained, it can be used
216directly ( TMultiLayerPerceptron::Evaluate() ) or exported to a
217standalone C++ code ( TMultiLayerPerceptron::Export() ).
218
219Finally, note that even if this implementation is inspired from the mlpfit code,
220the feature lists are not exactly matching:
221
222 - mlpfit hybrid learning method is not implemented
223 - output neurons can be normalized, this is not the case for mlpfit
224 - the neural net is exported in C++, FORTRAN or PYTHON
225 - the drawResult() method allows a fast check of the learning procedure
226
227In addition, the paw version of mlpfit had additional limitations on the number of
228neurons, hidden layers and inputs/outputs that does not apply to TMultiLayerPerceptron.
229*/
230
231
233#include "TSynapse.h"
234#include "TNeuron.h"
235#include "TClass.h"
236#include "TTree.h"
237#include "TEventList.h"
238#include "TRandom3.h"
239#include "TTimeStamp.h"
240#include "TRegexp.h"
241#include "TCanvas.h"
242#include "TH2.h"
243#include "TGraph.h"
244#include "TLegend.h"
245#include "TMultiGraph.h"
246#include "TDirectory.h"
247#include "TSystem.h"
248#include <iostream>
249#include <fstream>
250#include "TMath.h"
251#include "TTreeFormula.h"
252#include "TTreeFormulaManager.h"
253#include "TMarker.h"
254#include "TLine.h"
255#include "TText.h"
256#include "TObjString.h"
257#include <cstdlib>
258
260
261////////////////////////////////////////////////////////////////////////////////
262/// Default constructor
263
265{
266 if(!TClass::GetClass("TTreePlayer")) gSystem->Load("libTreePlayer");
267 fNetwork.SetOwner(true);
268 fFirstLayer.SetOwner(false);
269 fLastLayer.SetOwner(false);
270 fSynapses.SetOwner(true);
271 fData = nullptr;
272 fCurrentTree = -1;
274 fStructure = "";
275 fWeight = "1";
276 fTraining = nullptr;
277 fTrainingOwner = false;
278 fTest = nullptr;
279 fTestOwner = false;
280 fEventWeight = nullptr;
281 fManager = nullptr;
283 fEta = .1;
284 fEtaDecay = 1;
285 fDelta = 0;
286 fEpsilon = 0;
287 fTau = 3;
288 fLastAlpha = 0;
289 fReset = 50;
292 fextF = "";
293 fextD = "";
294}
295
296////////////////////////////////////////////////////////////////////////////////
297/// The network is described by a simple string:
298/// The input/output layers are defined by giving
299/// the branch names separated by comas.
300/// Hidden layers are just described by the number of neurons.
301/// The layers are separated by colons.
302///
303/// Ex: "x,y:10:5:f"
304///
305/// The output can be prepended by '@' if the variable has to be
306/// normalized.
307/// The output can be followed by '!' to use Softmax neurons for the
308/// output layer only.
309///
310/// Ex: "x,y:10:5:c1,c2,c3!"
311///
312/// Input and outputs are taken from the TTree given as second argument.
313/// training and test are the two TEventLists defining events
314/// to be used during the neural net training.
315/// Both the TTree and the TEventLists can be defined in the constructor,
316/// or later with the suited setter method.
317
319 TEventList * training,
320 TEventList * test,
322 const char* extF, const char* extD)
323{
324 if(!TClass::GetClass("TTreePlayer")) gSystem->Load("libTreePlayer");
325 fNetwork.SetOwner(true);
326 fFirstLayer.SetOwner(false);
327 fLastLayer.SetOwner(false);
328 fSynapses.SetOwner(true);
329 fStructure = layout;
330 fData = data;
331 fCurrentTree = -1;
333 fTraining = training;
334 fTrainingOwner = false;
335 fTest = test;
336 fTestOwner = false;
337 fWeight = "1";
338 fType = type;
340 fextF = extF;
341 fextD = extD;
342 fEventWeight = nullptr;
343 fManager = nullptr;
344 if (data) {
345 BuildNetwork();
346 AttachData();
347 }
349 fEta = .1;
350 fEpsilon = 0;
351 fDelta = 0;
352 fEtaDecay = 1;
353 fTau = 3;
354 fLastAlpha = 0;
355 fReset = 50;
356}
357
358////////////////////////////////////////////////////////////////////////////////
359/// The network is described by a simple string:
360/// The input/output layers are defined by giving
361/// the branch names separated by comas.
362/// Hidden layers are just described by the number of neurons.
363/// The layers are separated by colons.
364///
365/// Ex: "x,y:10:5:f"
366///
367/// The output can be prepended by '@' if the variable has to be
368/// normalized.
369/// The output can be followed by '!' to use Softmax neurons for the
370/// output layer only.
371///
372/// Ex: "x,y:10:5:c1,c2,c3!"
373///
374/// Input and outputs are taken from the TTree given as second argument.
375/// training and test are the two TEventLists defining events
376/// to be used during the neural net training.
377/// Both the TTree and the TEventLists can be defined in the constructor,
378/// or later with the suited setter method.
379
381 const char * weight, TTree * data,
382 TEventList * training,
383 TEventList * test,
385 const char* extF, const char* extD)
386{
387 if(!TClass::GetClass("TTreePlayer")) gSystem->Load("libTreePlayer");
388 fNetwork.SetOwner(true);
389 fFirstLayer.SetOwner(false);
390 fLastLayer.SetOwner(false);
391 fSynapses.SetOwner(true);
392 fStructure = layout;
393 fData = data;
394 fCurrentTree = -1;
396 fTraining = training;
397 fTrainingOwner = false;
398 fTest = test;
399 fTestOwner = false;
400 fWeight = weight;
401 fType = type;
403 fextF = extF;
404 fextD = extD;
405 fEventWeight = nullptr;
406 fManager = nullptr;
407 if (data) {
408 BuildNetwork();
409 AttachData();
410 }
412 fEta = .1;
413 fEtaDecay = 1;
414 fDelta = 0;
415 fEpsilon = 0;
416 fTau = 3;
417 fLastAlpha = 0;
418 fReset = 50;
419}
420
421////////////////////////////////////////////////////////////////////////////////
422/// The network is described by a simple string:
423/// The input/output layers are defined by giving
424/// the branch names separated by comas.
425/// Hidden layers are just described by the number of neurons.
426/// The layers are separated by colons.
427///
428/// Ex: "x,y:10:5:f"
429///
430/// The output can be prepended by '@' if the variable has to be
431/// normalized.
432/// The output can be followed by '!' to use Softmax neurons for the
433/// output layer only.
434///
435/// Ex: "x,y:10:5:c1,c2,c3!"
436///
437/// Input and outputs are taken from the TTree given as second argument.
438/// training and test are two cuts (see TTreeFormula) defining events
439/// to be used during the neural net training and testing.
440///
441/// Example: "Entry$%2", "(Entry$+1)%2".
442///
443/// Both the TTree and the cut can be defined in the constructor,
444/// or later with the suited setter method.
445
447 const char * training,
448 const char * test,
450 const char* extF, const char* extD)
451{
452 if(!TClass::GetClass("TTreePlayer")) gSystem->Load("libTreePlayer");
453 fNetwork.SetOwner(true);
454 fFirstLayer.SetOwner(false);
455 fLastLayer.SetOwner(false);
456 fSynapses.SetOwner(true);
457 fStructure = layout;
458 fData = data;
459 fCurrentTree = -1;
461 {
463 fTraining = new TEventList(Form("fTrainingList_%zu",(size_t)this));
464 }
465 fTrainingOwner = true;
466 {
468 fTest = new TEventList(Form("fTestList_%zu",(size_t)this));
469 }
470 fTestOwner = true;
471 fWeight = "1";
472 TString testcut = test;
473 if(testcut=="") testcut = Form("!(%s)",training);
474 fType = type;
476 fextF = extF;
477 fextD = extD;
478 fEventWeight = nullptr;
479 fManager = nullptr;
480 if (data) {
481 BuildNetwork();
482 data->Draw(Form(">>fTrainingList_%zu",(size_t)this),training,"goff");
483 data->Draw(Form(">>fTestList_%zu",(size_t)this),(const char *)testcut,"goff");
484 AttachData();
485 }
486 else {
487 Warning("TMultiLayerPerceptron::TMultiLayerPerceptron","Data not set. Cannot define datasets");
488 }
490 fEta = .1;
491 fEtaDecay = 1;
492 fDelta = 0;
493 fEpsilon = 0;
494 fTau = 3;
495 fLastAlpha = 0;
496 fReset = 50;
497}
498
499////////////////////////////////////////////////////////////////////////////////
500/// The network is described by a simple string:
501/// The input/output layers are defined by giving
502/// the branch names separated by comas.
503/// Hidden layers are just described by the number of neurons.
504/// The layers are separated by colons.
505///
506/// Ex: "x,y:10:5:f"
507///
508/// The output can be prepended by '@' if the variable has to be
509/// normalized.
510/// The output can be followed by '!' to use Softmax neurons for the
511/// output layer only.
512///
513/// Ex: "x,y:10:5:c1,c2,c3!"
514///
515/// Input and outputs are taken from the TTree given as second argument.
516/// training and test are two cuts (see TTreeFormula) defining events
517/// to be used during the neural net training and testing.
518///
519/// Example: "Entry$%2", "(Entry$+1)%2".
520///
521/// Both the TTree and the cut can be defined in the constructor,
522/// or later with the suited setter method.
523
525 const char * weight, TTree * data,
526 const char * training,
527 const char * test,
529 const char* extF, const char* extD)
530{
531 if(!TClass::GetClass("TTreePlayer")) gSystem->Load("libTreePlayer");
532 fNetwork.SetOwner(true);
533 fFirstLayer.SetOwner(false);
534 fLastLayer.SetOwner(false);
535 fSynapses.SetOwner(true);
536 fStructure = layout;
537 fData = data;
538 fCurrentTree = -1;
540 {
542 fTraining = new TEventList(Form("fTrainingList_%zu",(size_t)this));
543 }
544 fTrainingOwner = true;
545 {
547 fTest = new TEventList(Form("fTestList_%zu",(size_t)this));
548 }
549 fTestOwner = true;
550 fWeight = weight;
551 TString testcut = test;
552 if(testcut=="") testcut = Form("!(%s)",training);
553 fType = type;
555 fextF = extF;
556 fextD = extD;
557 fEventWeight = nullptr;
558 fManager = nullptr;
559 if (data) {
560 BuildNetwork();
561 data->Draw(Form(">>fTrainingList_%zu",(size_t)this),training,"goff");
562 data->Draw(Form(">>fTestList_%zu",(size_t)this),(const char *)testcut,"goff");
563 AttachData();
564 }
565 else {
566 Warning("TMultiLayerPerceptron::TMultiLayerPerceptron","Data not set. Cannot define datasets");
567 }
569 fEta = .1;
570 fEtaDecay = 1;
571 fDelta = 0;
572 fEpsilon = 0;
573 fTau = 3;
574 fLastAlpha = 0;
575 fReset = 50;
576}
577
578////////////////////////////////////////////////////////////////////////////////
579/// Destructor
580
582{
583 if(fTraining && fTrainingOwner) delete fTraining;
584 if(fTest && fTestOwner) delete fTest;
585}
586
587////////////////////////////////////////////////////////////////////////////////
588/// Set the data source
589
591{
592 if (fData) {
593 std::cerr << "Error: data already defined." << std::endl;
594 return;
595 }
596 fData = data;
597 if (data) {
598 BuildNetwork();
599 AttachData();
600 }
601}
602
603////////////////////////////////////////////////////////////////////////////////
604/// Set the event weight
605
607{
608 fWeight=branch;
609 if (fData) {
610 if (fEventWeight) {
612 delete fEventWeight;
613 }
614 fManager->Add((fEventWeight = new TTreeFormula("NNweight",fWeight.Data(),fData)));
615 }
616}
617
618////////////////////////////////////////////////////////////////////////////////
619/// Sets the Training dataset.
620/// Those events will be used for the minimization
621
623{
624 if(fTraining && fTrainingOwner) delete fTraining;
625 fTraining = train;
626 fTrainingOwner = false;
627}
628
629////////////////////////////////////////////////////////////////////////////////
630/// Sets the Test dataset.
631/// Those events will not be used for the minimization but for control
632
634{
635 if(fTest && fTestOwner) delete fTest;
636 fTest = test;
637 fTestOwner = false;
638}
639
640////////////////////////////////////////////////////////////////////////////////
641/// Sets the Training dataset.
642/// Those events will be used for the minimization.
643/// Note that the tree must be already defined.
644
646{
647 if(fTraining && fTrainingOwner) delete fTraining;
648 {
650 fTraining = new TEventList(Form("fTrainingList_%zu",(size_t)this));
651 }
652 fTrainingOwner = true;
653 if (fData) {
654 fData->Draw(Form(">>fTrainingList_%zu",(size_t)this),train,"goff");
655 }
656 else {
657 Warning("TMultiLayerPerceptron::TMultiLayerPerceptron","Data not set. Cannot define datasets");
658 }
659}
660
661////////////////////////////////////////////////////////////////////////////////
662/// Sets the Test dataset.
663/// Those events will not be used for the minimization but for control.
664/// Note that the tree must be already defined.
665
667{
668 if(fTest && fTestOwner) {delete fTest; fTest=nullptr;}
669 if(fTest) if(strncmp(fTest->GetName(),Form("fTestList_%zu",(size_t)this),10)) delete fTest;
670 {
672 fTest = new TEventList(Form("fTestList_%zu",(size_t)this));
673 }
674 fTestOwner = true;
675 if (fData) {
676 fData->Draw(Form(">>fTestList_%zu",(size_t)this),test,"goff");
677 }
678 else {
679 Warning("TMultiLayerPerceptron::TMultiLayerPerceptron","Data not set. Cannot define datasets");
680 }
681}
682
683////////////////////////////////////////////////////////////////////////////////
684/// Sets the learning method.
685/// Available methods are: kStochastic, kBatch,
686/// kSteepestDescent, kRibierePolak, kFletcherReeves and kBFGS.
687/// (look at the constructor for the complete description
688/// of learning methods and parameters)
689
691{
692 fLearningMethod = method;
693}
694
695////////////////////////////////////////////////////////////////////////////////
696/// Sets Eta - used in stochastic minimisation
697/// (look at the constructor for the complete description
698/// of learning methods and parameters)
699
701{
702 fEta = eta;
703}
704
705////////////////////////////////////////////////////////////////////////////////
706/// Sets Epsilon - used in stochastic minimisation
707/// (look at the constructor for the complete description
708/// of learning methods and parameters)
709
711{
712 fEpsilon = eps;
713}
714
715////////////////////////////////////////////////////////////////////////////////
716/// Sets Delta - used in stochastic minimisation
717/// (look at the constructor for the complete description
718/// of learning methods and parameters)
719
721{
722 fDelta = delta;
723}
724
725////////////////////////////////////////////////////////////////////////////////
726/// Sets EtaDecay - Eta *= EtaDecay at each epoch
727/// (look at the constructor for the complete description
728/// of learning methods and parameters)
729
731{
732 fEtaDecay = ed;
733}
734
735////////////////////////////////////////////////////////////////////////////////
736/// Sets Tau - used in line search
737/// (look at the constructor for the complete description
738/// of learning methods and parameters)
739
741{
742 fTau = tau;
743}
744
745////////////////////////////////////////////////////////////////////////////////
746/// Sets number of epochs between two resets of the
747/// search direction to the steepest descent.
748/// (look at the constructor for the complete description
749/// of learning methods and parameters)
750
752{
753 fReset = reset;
754}
755
756////////////////////////////////////////////////////////////////////////////////
757/// Load an entry into the network
758
760{
761 if (!fData) return;
762 fData->GetEntry(entry);
763 if (fData->GetTreeNumber() != fCurrentTree) {
764 ((TMultiLayerPerceptron*)this)->fCurrentTree = fData->GetTreeNumber();
765 fManager->Notify();
766 ((TMultiLayerPerceptron*)this)->fCurrentTreeWeight = fData->GetWeight();
767 }
769 for (Int_t i=0;i<nentries;i++) {
770 TNeuron *neuron = (TNeuron *)fNetwork.UncheckedAt(i);
771 neuron->SetNewEvent();
772 }
773}
774
775////////////////////////////////////////////////////////////////////////////////
776/// Train the network.
777/// nEpoch is the number of iterations.
778/// option can contain:
779/// - "text" (simple text output)
780/// - "graph" (evoluting graphical training curves)
781/// - "update=X" (step for the text/graph output update)
782/// - "+" will skip the randomisation and start from the previous values.
783/// - "current" (draw in the current canvas)
784/// - "minErrorTrain" (stop when NN error on the training sample gets below minE
785/// - "minErrorTest" (stop when NN error on the test sample gets below minE
786/// All combinations are available.
787
789{
790 Int_t i;
791 TString opt = option;
792 opt.ToLower();
793 // Decode options and prepare training.
794 Int_t verbosity = 0;
795 Bool_t newCanvas = true;
796 Bool_t minE_Train = false;
797 Bool_t minE_Test = false;
798 if (opt.Contains("text"))
799 verbosity += 1;
800 if (opt.Contains("graph"))
801 verbosity += 2;
802 Int_t displayStepping = 1;
803 if (opt.Contains("update=")) {
804 TRegexp reg("update=[0-9]*");
805 TString out = opt(reg);
806 displayStepping = atoi(out.Data() + 7);
807 }
808 if (opt.Contains("current"))
809 newCanvas = false;
810 if (opt.Contains("minerrortrain"))
811 minE_Train = true;
812 if (opt.Contains("minerrortest"))
813 minE_Test = true;
814 TVirtualPad *canvas = nullptr;
815 TMultiGraph *residual_plot = nullptr;
816 TGraph *train_residual_plot = nullptr;
817 TGraph *test_residual_plot = nullptr;
818 if ((!fData) || (!fTraining) || (!fTest)) {
819 Error("Train","Training/Test samples still not defined. Cannot train the neural network");
820 return;
821 }
822 Info("Train","Using %d train and %d test entries.",
823 fTraining->GetN(), fTest->GetN());
824 // Text and Graph outputs
825 if (verbosity % 2)
826 std::cout << "Training the Neural Network" << std::endl;
827 if (verbosity / 2) {
828 residual_plot = new TMultiGraph;
829 if(newCanvas)
830 canvas = new TCanvas("NNtraining", "Neural Net training");
831 else {
832 canvas = gPad;
833 if(!canvas) canvas = new TCanvas("NNtraining", "Neural Net training");
834 }
835 train_residual_plot = new TGraph(nEpoch);
836 test_residual_plot = new TGraph(nEpoch);
837 canvas->SetLeftMargin(0.14);
838 train_residual_plot->SetLineColor(4);
839 test_residual_plot->SetLineColor(2);
840 residual_plot->Add(train_residual_plot);
841 residual_plot->Add(test_residual_plot);
842 residual_plot->Draw("LA");
843 if (residual_plot->GetXaxis()) residual_plot->GetXaxis()->SetTitle("Epoch");
844 if (residual_plot->GetYaxis()) residual_plot->GetYaxis()->SetTitle("Error");
845 }
846 // If the option "+" is not set, one has to randomize the weights first
847 if (!opt.Contains("+"))
848 Randomize();
849 // Initialisation
850 fLastAlpha = 0;
852 Double_t *buffer = new Double_t[els];
853 Double_t *dir = new Double_t[els];
854 for (i = 0; i < els; i++)
855 buffer[i] = 0;
856 Int_t matrix_size = fLearningMethod==TMultiLayerPerceptron::kBFGS ? els : 1;
857 TMatrixD bfgsh(matrix_size, matrix_size);
858 TMatrixD gamma(matrix_size, 1);
859 TMatrixD delta(matrix_size, 1);
860 // Epoch loop. Here is the training itself.
861 Double_t training_E = 1e10;
862 Double_t test_E = 1e10;
863 for (Int_t iepoch = 0; (iepoch < nEpoch) && (!minE_Train || training_E>minE) && (!minE_Test || test_E>minE) ; iepoch++) {
864 switch (fLearningMethod) {
866 {
867 MLP_Stochastic(buffer);
868 break;
869 }
871 {
872 ComputeDEDw();
873 MLP_Batch(buffer);
874 break;
875 }
877 {
878 ComputeDEDw();
879 SteepestDir(dir);
880 if (LineSearch(dir, buffer))
881 MLP_Batch(buffer);
882 break;
883 }
885 {
886 ComputeDEDw();
887 if (!(iepoch % fReset)) {
888 SteepestDir(dir);
889 } else {
890 Double_t norm = 0;
891 Double_t onorm = 0;
892 for (i = 0; i < els; i++)
893 onorm += dir[i] * dir[i];
894 Double_t prod = 0;
895 Int_t idx = 0;
896 TNeuron *neuron = nullptr;
897 TSynapse *synapse = nullptr;
899 for (i=0;i<nentries;i++) {
900 neuron = (TNeuron *) fNetwork.UncheckedAt(i);
901 prod -= dir[idx++] * neuron->GetDEDw();
902 norm += neuron->GetDEDw() * neuron->GetDEDw();
903 }
905 for (i=0;i<nentries;i++) {
906 synapse = (TSynapse *) fSynapses.UncheckedAt(i);
907 prod -= dir[idx++] * synapse->GetDEDw();
908 norm += synapse->GetDEDw() * synapse->GetDEDw();
909 }
910 ConjugateGradientsDir(dir, (norm - prod) / onorm);
911 }
912 if (LineSearch(dir, buffer))
913 MLP_Batch(buffer);
914 break;
915 }
917 {
918 ComputeDEDw();
919 if (!(iepoch % fReset)) {
920 SteepestDir(dir);
921 } else {
922 Double_t norm = 0;
923 Double_t onorm = 0;
924 for (i = 0; i < els; i++)
925 onorm += dir[i] * dir[i];
926 TNeuron *neuron = nullptr;
927 TSynapse *synapse = nullptr;
929 for (i=0;i<nentries;i++) {
930 neuron = (TNeuron *) fNetwork.UncheckedAt(i);
931 norm += neuron->GetDEDw() * neuron->GetDEDw();
932 }
934 for (i=0;i<nentries;i++) {
935 synapse = (TSynapse *) fSynapses.UncheckedAt(i);
936 norm += synapse->GetDEDw() * synapse->GetDEDw();
937 }
938 ConjugateGradientsDir(dir, norm / onorm);
939 }
940 if (LineSearch(dir, buffer))
941 MLP_Batch(buffer);
942 break;
943 }
945 {
946 SetGammaDelta(gamma, delta, buffer);
947 if (!(iepoch % fReset)) {
948 SteepestDir(dir);
949 bfgsh.UnitMatrix();
950 } else {
951 if (GetBFGSH(bfgsh, gamma, delta)) {
952 SteepestDir(dir);
953 bfgsh.UnitMatrix();
954 } else {
955 BFGSDir(bfgsh, dir);
956 }
957 }
958 if (DerivDir(dir) > 0) {
959 SteepestDir(dir);
960 bfgsh.UnitMatrix();
961 }
962 if (LineSearch(dir, buffer)) {
963 bfgsh.UnitMatrix();
964 SteepestDir(dir);
965 if (LineSearch(dir, buffer)) {
966 Error("TMultiLayerPerceptron::Train()","Line search fail");
967 iepoch = nEpoch;
968 }
969 }
970 break;
971 }
972 }
973 // Security: would the learning lead to non real numbers,
974 // the learning should stop now.
976 Error("TMultiLayerPerceptron::Train()","Stop.");
977 iepoch = nEpoch;
978 }
979 // Process other ROOT events. Time penalty is less than
980 // 1/1000 sec/evt on a mobile AMD Athlon(tm) XP 1500+
984 // Intermediate graph and text output
985 if ((verbosity % 2) && ((!(iepoch % displayStepping)) || (iepoch == nEpoch - 1))) {
986 std::cout << "Epoch: " << iepoch
987 << " learn=" << training_E
988 << " test=" << test_E
989 << std::endl;
990 }
991 if (verbosity / 2) {
992 train_residual_plot->SetPoint(iepoch, iepoch,training_E);
993 test_residual_plot->SetPoint(iepoch, iepoch,test_E);
994 if (!iepoch) {
995 Double_t trp = train_residual_plot->GetY()[iepoch];
996 Double_t tep = test_residual_plot->GetY()[iepoch];
997 for (i = 1; i < nEpoch; i++) {
998 train_residual_plot->SetPoint(i, i, trp);
999 test_residual_plot->SetPoint(i, i, tep);
1000 }
1001 }
1002 if ((!(iepoch % displayStepping)) || (iepoch == nEpoch - 1)) {
1003 if (residual_plot->GetYaxis()) {
1004 residual_plot->GetYaxis()->UnZoom();
1005 residual_plot->GetYaxis()->SetTitleOffset(1.4);
1006 residual_plot->GetYaxis()->SetDecimals();
1007 }
1008 canvas->Modified();
1009 canvas->Update();
1010 }
1011 }
1012 }
1013 // Cleaning
1014 delete [] buffer;
1015 delete [] dir;
1016 // Final Text and Graph outputs
1017 if (verbosity % 2)
1018 std::cout << "Training done." << std::endl;
1019 if (verbosity / 2) {
1020 TLegend *legend = new TLegend(.75, .80, .95, .95);
1021 legend->AddEntry(residual_plot->GetListOfGraphs()->At(0),
1022 "Training sample", "L");
1023 legend->AddEntry(residual_plot->GetListOfGraphs()->At(1),
1024 "Test sample", "L");
1025 legend->Draw();
1026 }
1027}
1028
1029////////////////////////////////////////////////////////////////////////////////
1030/// Computes the output for a given event.
1031/// Look at the output neuron designed by index.
1032
1034{
1035 GetEntry(event);
1036 TNeuron *out = (TNeuron *) (fLastLayer.At(index));
1037 if (out)
1038 return out->GetValue();
1039 else
1040 return 0;
1041}
1042
1043////////////////////////////////////////////////////////////////////////////////
1044/// Error on the output for a given event
1045
1047{
1048 GetEntry(event);
1049 Double_t error = 0;
1050 // look at 1st output neuron to determine type and error function
1051 Int_t nEntries = fLastLayer.GetEntriesFast();
1052 if (nEntries == 0) return 0.0;
1053 switch (fOutType) {
1054 case (TNeuron::kSigmoid):
1055 error = GetCrossEntropyBinary();
1056 break;
1057 case (TNeuron::kSoftmax):
1058 error = GetCrossEntropy();
1059 break;
1060 case (TNeuron::kLinear):
1061 error = GetSumSquareError();
1062 break;
1063 default:
1064 // default to sum-of-squares error
1065 error = GetSumSquareError();
1066 }
1067 error *= fEventWeight->EvalInstance();
1068 error *= fCurrentTreeWeight;
1069 return error;
1070}
1071
1072////////////////////////////////////////////////////////////////////////////////
1073/// Error on the whole dataset
1074
1076{
1077 TEventList *list =
1079 Double_t error = 0;
1080 Int_t i;
1081 if (list) {
1082 Int_t nEvents = list->GetN();
1083 for (i = 0; i < nEvents; i++) {
1084 error += GetError(list->GetEntry(i));
1085 }
1086 } else if (fData) {
1087 Int_t nEvents = (Int_t) fData->GetEntries();
1088 for (i = 0; i < nEvents; i++) {
1089 error += GetError(i);
1090 }
1091 }
1092 return error;
1093}
1094
1095////////////////////////////////////////////////////////////////////////////////
1096/// Error on the output for a given event
1097
1099{
1100 Double_t error = 0;
1101 for (Int_t i = 0; i < fLastLayer.GetEntriesFast(); i++) {
1102 TNeuron *neuron = (TNeuron *) fLastLayer[i];
1103 error += neuron->GetError() * neuron->GetError();
1104 }
1105 return (error / 2.);
1106}
1107
1108////////////////////////////////////////////////////////////////////////////////
1109/// Cross entropy error for sigmoid output neurons, for a given event
1110
1112{
1113 Double_t error = 0;
1114 for (Int_t i = 0; i < fLastLayer.GetEntriesFast(); i++) {
1115 TNeuron *neuron = (TNeuron *) fLastLayer[i];
1116 Double_t output = neuron->GetValue(); // sigmoid output and target
1117 Double_t target = neuron->GetTarget(); // values lie in [0,1]
1118 if (target < DBL_EPSILON) {
1119 if (output == 1.0)
1120 error = DBL_MAX;
1121 else
1122 error -= TMath::Log(1 - output);
1123 } else
1124 if ((1 - target) < DBL_EPSILON) {
1125 if (output == 0.0)
1126 error = DBL_MAX;
1127 else
1128 error -= TMath::Log(output);
1129 } else {
1130 if (output == 0.0 || output == 1.0)
1131 error = DBL_MAX;
1132 else
1133 error -= target * TMath::Log(output / target) + (1-target) * TMath::Log((1 - output)/(1 - target));
1134 }
1135 }
1136 return error;
1137}
1138
1139////////////////////////////////////////////////////////////////////////////////
1140/// Cross entropy error for a softmax output neuron, for a given event
1141
1143{
1144 Double_t error = 0;
1145 for (Int_t i = 0; i < fLastLayer.GetEntriesFast(); i++) {
1146 TNeuron *neuron = (TNeuron *) fLastLayer[i];
1147 Double_t output = neuron->GetValue(); // softmax output and target
1148 Double_t target = neuron->GetTarget(); // values lie in [0,1]
1149 if (target > DBL_EPSILON) { // (target == 0) => dE = 0
1150 if (output == 0.0)
1151 error = DBL_MAX;
1152 else
1153 error -= target * TMath::Log(output / target);
1154 }
1155 }
1156 return error;
1157}
1158
1159////////////////////////////////////////////////////////////////////////////////
1160/// Compute the DEDw = sum on all training events of dedw for each weight
1161/// normalized by the number of events.
1162
1164{
1165 Int_t i,j;
1167 TSynapse *synapse;
1168 for (i=0;i<nentries;i++) {
1169 synapse = (TSynapse *) fSynapses.UncheckedAt(i);
1170 synapse->SetDEDw(0.);
1171 }
1172 TNeuron *neuron;
1174 for (i=0;i<nentries;i++) {
1175 neuron = (TNeuron *) fNetwork.UncheckedAt(i);
1176 neuron->SetDEDw(0.);
1177 }
1178 Double_t eventWeight = 1.;
1179 if (fTraining) {
1180 Int_t nEvents = fTraining->GetN();
1181 for (i = 0; i < nEvents; i++) {
1183 eventWeight = fEventWeight->EvalInstance();
1184 eventWeight *= fCurrentTreeWeight;
1186 for (j=0;j<nentries;j++) {
1187 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
1188 synapse->SetDEDw(synapse->GetDEDw() + (synapse->GetDeDw()*eventWeight));
1189 }
1191 for (j=0;j<nentries;j++) {
1192 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
1193 neuron->SetDEDw(neuron->GetDEDw() + (neuron->GetDeDw()*eventWeight));
1194 }
1195 }
1197 for (j=0;j<nentries;j++) {
1198 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
1199 synapse->SetDEDw(synapse->GetDEDw() / (Double_t) nEvents);
1200 }
1202 for (j=0;j<nentries;j++) {
1203 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
1204 neuron->SetDEDw(neuron->GetDEDw() / (Double_t) nEvents);
1205 }
1206 } else if (fData) {
1207 Int_t nEvents = (Int_t) fData->GetEntries();
1208 for (i = 0; i < nEvents; i++) {
1209 GetEntry(i);
1210 eventWeight = fEventWeight->EvalInstance();
1211 eventWeight *= fCurrentTreeWeight;
1213 for (j=0;j<nentries;j++) {
1214 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
1215 synapse->SetDEDw(synapse->GetDEDw() + (synapse->GetDeDw()*eventWeight));
1216 }
1218 for (j=0;j<nentries;j++) {
1219 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
1220 neuron->SetDEDw(neuron->GetDEDw() + (neuron->GetDeDw()*eventWeight));
1221 }
1222 }
1224 for (j=0;j<nentries;j++) {
1225 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
1226 synapse->SetDEDw(synapse->GetDEDw() / (Double_t) nEvents);
1227 }
1229 for (j=0;j<nentries;j++) {
1230 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
1231 neuron->SetDEDw(neuron->GetDEDw() / (Double_t) nEvents);
1232 }
1233 }
1234}
1235
1236////////////////////////////////////////////////////////////////////////////////
1237/// Randomize the weights
1238
1240{
1242 Int_t j;
1243 TSynapse *synapse;
1244 TNeuron *neuron;
1245 TTimeStamp ts;
1246 TRandom3 gen(ts.GetSec());
1247 for (j=0;j<nentries;j++) {
1248 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
1249 synapse->SetWeight(gen.Rndm() - 0.5);
1250 }
1252 for (j=0;j<nentries;j++) {
1253 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
1254 neuron->SetWeight(gen.Rndm() - 0.5);
1255 }
1256}
1257
1258////////////////////////////////////////////////////////////////////////////////
1259/// Connects the TTree to Neurons in input and output
1260/// layers. The formulas associated to each neuron are created
1261/// and reported to the network formula manager.
1262/// By default, the branch is not normalised since this would degrade
1263/// performance for classification jobs.
1264/// Normalisation can be requested by putting '@' in front of the formula.
1265
1267{
1268 Int_t j = 0;
1269 TNeuron *neuron = nullptr;
1270 Bool_t normalize = false;
1272
1273 // Set the size of the internal array of parameters of the formula
1274 Int_t maxop, maxpar, maxconst;
1275 ROOT::v5::TFormula::GetMaxima(maxop, maxpar, maxconst);
1277
1278 //first layer
1279 const TString input = TString(fStructure(0, fStructure.First(':')));
1280 const TObjArray *inpL = input.Tokenize(", ");
1282 // make sure nentries == entries in inpL
1283 R__ASSERT(nentries == inpL->GetLast()+1);
1284 for (j=0;j<nentries;j++) {
1285 normalize = false;
1286 const TString brName = ((TObjString *)inpL->At(j))->GetString();
1287 neuron = (TNeuron *) fFirstLayer.UncheckedAt(j);
1288 if (brName[0]=='@')
1289 normalize = true;
1290 fManager->Add(neuron->UseBranch(fData,brName.Data() + (normalize?1:0)));
1291 if(!normalize) neuron->SetNormalisation(0., 1.);
1292 }
1293 delete inpL;
1294
1295 // last layer
1297 fStructure(fStructure.Last(':') + 1,
1298 fStructure.Length() - fStructure.Last(':')));
1299 const TObjArray *outL = output.Tokenize(", ");
1301 // make sure nentries == entries in outL
1302 R__ASSERT(nentries == outL->GetLast()+1);
1303 for (j=0;j<nentries;j++) {
1304 normalize = false;
1305 const TString brName = ((TObjString *)outL->At(j))->GetString();
1306 neuron = (TNeuron *) fLastLayer.UncheckedAt(j);
1307 if (brName[0]=='@')
1308 normalize = true;
1309 fManager->Add(neuron->UseBranch(fData,brName.Data() + (normalize?1:0)));
1310 if(!normalize) neuron->SetNormalisation(0., 1.);
1311 }
1312 delete outL;
1313
1314 fManager->Add((fEventWeight = new TTreeFormula("NNweight",fWeight.Data(),fData)));
1315 //fManager->Sync();
1316
1317 // Set the old values
1318 ROOT::v5::TFormula::SetMaxima(maxop, maxpar, maxconst);
1319}
1320
1321////////////////////////////////////////////////////////////////////////////////
1322/// Expand the structure of the first layer
1323
1325{
1327 const TObjArray *inpL = input.Tokenize(", ");
1328 Int_t nneurons = inpL->GetLast()+1;
1329
1330 TString hiddenAndOutput = TString(
1331 fStructure(fStructure.First(':') + 1,
1332 fStructure.Length() - fStructure.First(':')));
1333 TString newInput;
1334 Int_t i = 0;
1335 // loop on input neurons
1336 for (i = 0; i<nneurons; i++) {
1337 const TString name = ((TObjString *)inpL->At(i))->GetString();
1338 TTreeFormula f("sizeTestFormula",name,fData);
1339 // Variable size arrays are unreliable
1340 if(f.GetMultiplicity()==1 && f.GetNdata()>1) {
1341 Warning("TMultiLayerPerceptron::ExpandStructure()","Variable size arrays cannot be used to build implicitly an input layer. The index 0 will be assumed.");
1342 }
1343 // Check if we are coping with an array... then expand
1344 // The array operator used is {}. It is detected in TNeuron, and
1345 // passed directly as instance index of the TTreeFormula,
1346 // so that complex compounds made of arrays can be used without
1347 // parsing the details.
1348 else if(f.GetNdata()>1) {
1349 for(Int_t j=0; j<f.GetNdata(); j++) {
1350 if(i||j) newInput += ",";
1351 newInput += name;
1352 newInput += "{";
1353 newInput += j;
1354 newInput += "}";
1355 }
1356 continue;
1357 }
1358 if(i) newInput += ",";
1359 newInput += name;
1360 }
1361 delete inpL;
1362
1363 // Save the result
1364 fStructure = newInput + ":" + hiddenAndOutput;
1365}
1366
1367////////////////////////////////////////////////////////////////////////////////
1368/// Instantiates the network from the description
1369
1371{
1374 TString hidden = TString(
1375 fStructure(fStructure.First(':') + 1,
1376 fStructure.Last(':') - fStructure.First(':') - 1));
1378 fStructure(fStructure.Last(':') + 1,
1379 fStructure.Length() - fStructure.Last(':')));
1380 Int_t bll = atoi(TString(
1381 hidden(hidden.Last(':') + 1,
1382 hidden.Length() - (hidden.Last(':') + 1))).Data());
1383 if (input.Length() == 0) {
1384 Error("BuildNetwork()","malformed structure. No input layer.");
1385 return;
1386 }
1387 if (output.Length() == 0) {
1388 Error("BuildNetwork()","malformed structure. No output layer.");
1389 return;
1390 }
1392 BuildHiddenLayers(hidden);
1393 BuildLastLayer(output, bll);
1394}
1395
1396////////////////////////////////////////////////////////////////////////////////
1397/// Instantiates the neurons in input
1398/// Inputs are normalised and the type is set to kOff
1399/// (simple forward of the formula value)
1400
1402{
1403 const TObjArray *inpL = input.Tokenize(", ");
1404 const Int_t nneurons =inpL->GetLast()+1;
1405 TNeuron *neuron = nullptr;
1406 Int_t i = 0;
1407 for (i = 0; i<nneurons; i++) {
1408 const TString name = ((TObjString *)inpL->At(i))->GetString();
1409 neuron = new TNeuron(TNeuron::kOff, name);
1410 fFirstLayer.AddLast(neuron);
1411 fNetwork.AddLast(neuron);
1412 }
1413 delete inpL;
1414}
1415
1416////////////////////////////////////////////////////////////////////////////////
1417/// Builds hidden layers.
1418
1420{
1421 Int_t beg = 0;
1422 Int_t end = hidden.Index(":", beg + 1);
1423 Int_t prevStart = 0;
1424 Int_t prevStop = fNetwork.GetEntriesFast();
1425 Int_t layer = 1;
1426 while (end != -1) {
1427 BuildOneHiddenLayer(hidden(beg, end - beg), layer, prevStart, prevStop, false);
1428 beg = end + 1;
1429 end = hidden.Index(":", beg + 1);
1430 }
1431
1432 BuildOneHiddenLayer(hidden(beg, hidden.Length() - beg), layer, prevStart, prevStop, true);
1433}
1434
1435////////////////////////////////////////////////////////////////////////////////
1436/// Builds a hidden layer, updates the number of layers.
1437
1439 Int_t& prevStart, Int_t& prevStop,
1440 Bool_t lastLayer)
1441{
1442 TNeuron *neuron = nullptr;
1443 TSynapse *synapse = nullptr;
1444 TString name;
1445 if (!sNumNodes.IsAlnum() || sNumNodes.IsAlpha()) {
1446 Error("BuildOneHiddenLayer",
1447 "The specification '%s' for hidden layer %d must contain only numbers!",
1448 sNumNodes.Data(), layer - 1);
1449 } else {
1450 Int_t num = atoi(sNumNodes.Data());
1451 for (Int_t i = 0; i < num; i++) {
1452 name.Form("HiddenL%d:N%d",layer,i);
1453 neuron = new TNeuron(fType, name, "", (const char*)fextF, (const char*)fextD);
1454 fNetwork.AddLast(neuron);
1455 for (Int_t j = prevStart; j < prevStop; j++) {
1456 synapse = new TSynapse((TNeuron *) fNetwork[j], neuron);
1457 fSynapses.AddLast(synapse);
1458 }
1459 }
1460
1461 if (!lastLayer) {
1462 // tell each neuron which ones are in its own layer (for Softmax)
1463 Int_t nEntries = fNetwork.GetEntriesFast();
1464 for (Int_t i = prevStop; i < nEntries; i++) {
1465 neuron = (TNeuron *) fNetwork[i];
1466 for (Int_t j = prevStop; j < nEntries; j++)
1467 neuron->AddInLayer((TNeuron *) fNetwork[j]);
1468 }
1469 }
1470
1471 prevStart = prevStop;
1472 prevStop = fNetwork.GetEntriesFast();
1473 layer++;
1474 }
1475}
1476
1477////////////////////////////////////////////////////////////////////////////////
1478/// Builds the output layer
1479/// Neurons are linear combinations of input, by default.
1480/// If the structure ends with "!", neurons are set up for classification,
1481/// ie. with a sigmoid (1 neuron) or softmax (more neurons) activation function.
1482
1484{
1485 Int_t nneurons = output.CountChar(',')+1;
1486 if (fStructure.EndsWith("!")) {
1487 fStructure = TString(fStructure(0, fStructure.Length() - 1)); // remove "!"
1488 if (nneurons == 1)
1490 else
1492 }
1493 Int_t prevStop = fNetwork.GetEntriesFast();
1494 Int_t prevStart = prevStop - prev;
1495 Ssiz_t pos = 0;
1496 TNeuron *neuron;
1497 TSynapse *synapse;
1498 TString name;
1499 Int_t i,j;
1500 for (i = 0; i<nneurons; i++) {
1501 Ssiz_t nextpos=output.Index(",",pos);
1502 if (nextpos!=kNPOS)
1503 name=output(pos,nextpos-pos);
1504 else name=output(pos,output.Length());
1505 pos+=nextpos+1;
1506 neuron = new TNeuron(fOutType, name);
1507 for (j = prevStart; j < prevStop; j++) {
1508 synapse = new TSynapse((TNeuron *) fNetwork[j], neuron);
1509 fSynapses.AddLast(synapse);
1510 }
1511 fLastLayer.AddLast(neuron);
1512 fNetwork.AddLast(neuron);
1513 }
1514 // tell each neuron which ones are in its own layer (for Softmax)
1515 Int_t nEntries = fNetwork.GetEntriesFast();
1516 for (i = prevStop; i < nEntries; i++) {
1517 neuron = (TNeuron *) fNetwork[i];
1518 for (j = prevStop; j < nEntries; j++)
1519 neuron->AddInLayer((TNeuron *) fNetwork[j]);
1520 }
1521
1522}
1523
1524////////////////////////////////////////////////////////////////////////////////
1525/// Draws the neural net output
1526/// It produces an histogram with the output for the two datasets.
1527/// Index is the number of the desired output neuron.
1528/// "option" can contain:
1529/// - test or train to select a dataset
1530/// - comp to produce a X-Y comparison plot
1531/// - nocanv to not create a new TCanvas for the plot
1532
1534{
1535 TString opt = option;
1536 opt.ToLower();
1537 TNeuron *out = (TNeuron *) (fLastLayer.At(index));
1538 if (!out) {
1539 Error("DrawResult()","no such output.");
1540 return;
1541 }
1542 //TCanvas *canvas = new TCanvas("NNresult", "Neural Net output");
1543 if (!opt.Contains("nocanv"))
1544 new TCanvas("NNresult", "Neural Net output");
1545 const Double_t *norm = out->GetNormalisation();
1546 TEventList *events = nullptr;
1547 TString setname;
1548 Int_t i;
1549 if (opt.Contains("train")) {
1550 events = fTraining;
1551 setname = Form("train%d",index);
1552 } else if (opt.Contains("test")) {
1553 events = fTest;
1554 setname = Form("test%d",index);
1555 }
1556 if ((!fData) || (!events)) {
1557 Error("DrawResult()","no dataset.");
1558 return;
1559 }
1560 if (opt.Contains("comp")) {
1561 //comparison plot
1562 TString title = "Neural Net Output control. ";
1563 title += setname;
1564 setname = "MLP_" + setname + "_comp";
1565 TH2D *hist = ((TH2D *) gDirectory->Get(setname.Data()));
1566 if (!hist)
1567 hist = new TH2D(setname.Data(), title.Data(), 50, -1, 1, 50, -1, 1);
1568 hist->Reset();
1569 Int_t nEvents = events->GetN();
1570 for (i = 0; i < nEvents; i++) {
1571 GetEntry(events->GetEntry(i));
1572 hist->Fill(out->GetValue(), (out->GetBranch() - norm[1]) / norm[0]);
1573 }
1574 hist->Draw();
1575 } else {
1576 //output plot
1577 TString title = "Neural Net Output. ";
1578 title += setname;
1579 setname = "MLP_" + setname;
1580 TH1D *hist = ((TH1D *) gDirectory->Get(setname.Data()));
1581 if (!hist)
1582 hist = new TH1D(setname, title, 50, 1, -1);
1583 hist->Reset();
1584 Int_t nEvents = events->GetN();
1585 for (i = 0; i < nEvents; i++)
1586 hist->Fill(Result(events->GetEntry(i), index));
1587 hist->Draw();
1588 if (opt.Contains("train") && opt.Contains("test")) {
1589 events = fTraining;
1590 setname = "train";
1591 hist = ((TH1D *) gDirectory->Get("MLP_test"));
1592 if (!hist)
1593 hist = new TH1D(setname, title, 50, 1, -1);
1594 hist->Reset();
1595 nEvents = events->GetN();
1596 for (i = 0; i < nEvents; i++)
1597 hist->Fill(Result(events->GetEntry(i), index));
1598 hist->Draw("same");
1599 }
1600 }
1601}
1602
1603////////////////////////////////////////////////////////////////////////////////
1604/// Dumps the weights to a text file.
1605/// Set filename to "-" (default) to dump to the standard output
1606
1608{
1609 TString filen = filename;
1610 std::ostream * output;
1611 if (filen == "") {
1612 Error("TMultiLayerPerceptron::DumpWeights()","Invalid file name");
1613 return kFALSE;
1614 }
1615 if (filen == "-")
1616 output = &std::cout;
1617 else
1618 output = new std::ofstream(filen.Data());
1619 TNeuron *neuron = nullptr;
1620 *output << "#input normalization" << std::endl;
1622 Int_t j=0;
1623 for (j=0;j<nentries;j++) {
1624 neuron = (TNeuron *) fFirstLayer.UncheckedAt(j);
1625 *output << neuron->GetNormalisation()[0] << " "
1626 << neuron->GetNormalisation()[1] << std::endl;
1627 }
1628 *output << "#output normalization" << std::endl;
1630 for (j=0;j<nentries;j++) {
1631 neuron = (TNeuron *) fLastLayer.UncheckedAt(j);
1632 *output << neuron->GetNormalisation()[0] << " "
1633 << neuron->GetNormalisation()[1] << std::endl;
1634 }
1635 *output << "#neurons weights" << std::endl;
1637 while ((neuron = (TNeuron *) it->Next()))
1638 *output << neuron->GetWeight() << std::endl;
1639 delete it;
1641 TSynapse *synapse = nullptr;
1642 *output << "#synapses weights" << std::endl;
1643 while ((synapse = (TSynapse *) it->Next()))
1644 *output << synapse->GetWeight() << std::endl;
1645 delete it;
1646 if (filen != "-") {
1647 ((std::ofstream *) output)->close();
1648 delete output;
1649 }
1650 return kTRUE;
1651}
1652
1653////////////////////////////////////////////////////////////////////////////////
1654/// Loads the weights from a text file conforming to the format
1655/// defined by DumpWeights.
1656
1658{
1659 TString filen = filename;
1660 Double_t w;
1661 if (filen == "") {
1662 Error("TMultiLayerPerceptron::LoadWeights()","Invalid file name");
1663 return kFALSE;
1664 }
1665 char *buff = new char[100];
1666 std::ifstream input(filen.Data());
1667 // input normalzation
1668 input.getline(buff, 100);
1670 Float_t n1,n2;
1671 TNeuron *neuron = nullptr;
1672 while ((neuron = (TNeuron *) it->Next())) {
1673 input >> n1 >> n2;
1674 neuron->SetNormalisation(n2,n1);
1675 }
1676 input.getline(buff, 100);
1677 // output normalization
1678 input.getline(buff, 100);
1679 delete it;
1681 while ((neuron = (TNeuron *) it->Next())) {
1682 input >> n1 >> n2;
1683 neuron->SetNormalisation(n2,n1);
1684 }
1685 input.getline(buff, 100);
1686 // neuron weights
1687 input.getline(buff, 100);
1688 delete it;
1690 while ((neuron = (TNeuron *) it->Next())) {
1691 input >> w;
1692 neuron->SetWeight(w);
1693 }
1694 delete it;
1695 input.getline(buff, 100);
1696 // synapse weights
1697 input.getline(buff, 100);
1699 TSynapse *synapse = nullptr;
1700 while ((synapse = (TSynapse *) it->Next())) {
1701 input >> w;
1702 synapse->SetWeight(w);
1703 }
1704 delete it;
1705 delete[] buff;
1706 return kTRUE;
1707}
1708
1709////////////////////////////////////////////////////////////////////////////////
1710/// Returns the Neural Net for a given set of input parameters
1711/// #%parameters must equal #%input neurons
1712
1714{
1716 TNeuron *neuron;
1717 while ((neuron = (TNeuron *) it->Next()))
1718 neuron->SetNewEvent();
1719 delete it;
1721 Int_t i=0;
1722 while ((neuron = (TNeuron *) it->Next()))
1723 neuron->ForceExternalValue(params[i++]);
1724 delete it;
1725 TNeuron *out = (TNeuron *) (fLastLayer.At(index));
1726 if (out)
1727 return out->GetValue();
1728 else
1729 return 0;
1730}
1731
1732////////////////////////////////////////////////////////////////////////////////
1733/// Exports the NN as a function for any non-ROOT-dependant code
1734/// Supported languages are: only C++ , FORTRAN and Python (yet)
1735/// This feature is also useful if you want to plot the NN as
1736/// a function (TF1 or TF2).
1737
1739{
1740 TString lg = language;
1741 lg.ToUpper();
1742 Int_t i;
1744 Warning("TMultiLayerPerceptron::Export","Request to export a network using an external function");
1745 }
1746 if (lg == "C++") {
1747 TString basefilename = filename;
1748 Int_t slash = basefilename.Last('/')+1;
1749 if (slash) basefilename = TString(basefilename(slash, basefilename.Length()-slash));
1750
1751 TString classname = basefilename;
1752 TString header = filename;
1753 header += ".h";
1754 TString source = filename;
1755 source += ".cxx";
1756 std::ofstream headerfile(header);
1757 std::ofstream sourcefile(source);
1758 headerfile << "#ifndef " << basefilename << "_h" << std::endl;
1759 headerfile << "#define " << basefilename << "_h" << std::endl << std::endl;
1760 headerfile << "class " << classname << " { " << std::endl;
1761 headerfile << "public:" << std::endl;
1762 headerfile << " " << classname << "() {}" << std::endl;
1763 headerfile << " ~" << classname << "() {}" << std::endl;
1764 sourcefile << "#include \"" << header << "\"" << std::endl;
1765 sourcefile << "#include <cmath>" << std::endl << std::endl;
1766 headerfile << " double Value(int index";
1767 sourcefile << "double " << classname << "::Value(int index";
1768 for (i = 0; i < fFirstLayer.GetEntriesFast(); i++) {
1769 headerfile << ",double in" << i;
1770 sourcefile << ",double in" << i;
1771 }
1772 headerfile << ");" << std::endl;
1773 sourcefile << ") {" << std::endl;
1774 for (i = 0; i < fFirstLayer.GetEntriesFast(); i++)
1775 sourcefile << " input" << i << " = (in" << i << " - "
1776 << ((TNeuron *) fFirstLayer[i])->GetNormalisation()[1] << ")/"
1777 << ((TNeuron *) fFirstLayer[i])->GetNormalisation()[0] << ";"
1778 << std::endl;
1779 sourcefile << " switch(index) {" << std::endl;
1780 TNeuron *neuron;
1782 Int_t idx = 0;
1783 while ((neuron = (TNeuron *) it->Next()))
1784 sourcefile << " case " << idx++ << ":" << std::endl
1785 << " return neuron" << neuron << "();" << std::endl;
1786 sourcefile << " default:" << std::endl
1787 << " return 0.;" << std::endl << " }"
1788 << std::endl;
1789 sourcefile << "}" << std::endl << std::endl;
1790 headerfile << " double Value(int index, double* input);" << std::endl;
1791 sourcefile << "double " << classname << "::Value(int index, double* input) {" << std::endl;
1792 for (i = 0; i < fFirstLayer.GetEntriesFast(); i++)
1793 sourcefile << " input" << i << " = (input[" << i << "] - "
1794 << ((TNeuron *) fFirstLayer[i])->GetNormalisation()[1] << ")/"
1795 << ((TNeuron *) fFirstLayer[i])->GetNormalisation()[0] << ";"
1796 << std::endl;
1797 sourcefile << " switch(index) {" << std::endl;
1798 delete it;
1800 idx = 0;
1801 while ((neuron = (TNeuron *) it->Next()))
1802 sourcefile << " case " << idx++ << ":" << std::endl
1803 << " return neuron" << neuron << "();" << std::endl;
1804 sourcefile << " default:" << std::endl
1805 << " return 0.;" << std::endl << " }"
1806 << std::endl;
1807 sourcefile << "}" << std::endl << std::endl;
1808 headerfile << "private:" << std::endl;
1809 for (i = 0; i < fFirstLayer.GetEntriesFast(); i++)
1810 headerfile << " double input" << i << ";" << std::endl;
1811 delete it;
1813 idx = 0;
1814 while ((neuron = (TNeuron *) it->Next())) {
1815 if (!neuron->GetPre(0)) {
1816 headerfile << " double neuron" << neuron << "();" << std::endl;
1817 sourcefile << "double " << classname << "::neuron" << neuron
1818 << "() {" << std::endl;
1819 sourcefile << " return input" << idx++ << ";" << std::endl;
1820 sourcefile << "}" << std::endl << std::endl;
1821 } else {
1822 headerfile << " double input" << neuron << "();" << std::endl;
1823 sourcefile << "double " << classname << "::input" << neuron
1824 << "() {" << std::endl;
1825 sourcefile << " double input = " << neuron->GetWeight()
1826 << ";" << std::endl;
1827 TSynapse *syn = nullptr;
1828 Int_t n = 0;
1829 while ((syn = neuron->GetPre(n++))) {
1830 sourcefile << " input += synapse" << syn << "();" << std::endl;
1831 }
1832 sourcefile << " return input;" << std::endl;
1833 sourcefile << "}" << std::endl << std::endl;
1834
1835 headerfile << " double neuron" << neuron << "();" << std::endl;
1836 sourcefile << "double " << classname << "::neuron" << neuron << "() {" << std::endl;
1837 sourcefile << " double input = input" << neuron << "();" << std::endl;
1838 switch(neuron->GetType()) {
1839 case (TNeuron::kSigmoid):
1840 {
1841 sourcefile << " return ((input < -709. ? 0. : (1/(1+exp(-input)))) * ";
1842 break;
1843 }
1844 case (TNeuron::kLinear):
1845 {
1846 sourcefile << " return (input * ";
1847 break;
1848 }
1849 case (TNeuron::kTanh):
1850 {
1851 sourcefile << " return (tanh(input) * ";
1852 break;
1853 }
1854 case (TNeuron::kGauss):
1855 {
1856 sourcefile << " return (exp(-input*input) * ";
1857 break;
1858 }
1859 case (TNeuron::kSoftmax):
1860 {
1861 sourcefile << " return (exp(input) / (";
1862 Int_t nn = 0;
1863 TNeuron* side = neuron->GetInLayer(nn++);
1864 sourcefile << "exp(input" << side << "())";
1865 while ((side = neuron->GetInLayer(nn++)))
1866 sourcefile << " + exp(input" << side << "())";
1867 sourcefile << ") * ";
1868 break;
1869 }
1870 default:
1871 {
1872 sourcefile << " return (0.0 * ";
1873 }
1874 }
1875 sourcefile << neuron->GetNormalisation()[0] << ")+" ;
1876 sourcefile << neuron->GetNormalisation()[1] << ";" << std::endl;
1877 sourcefile << "}" << std::endl << std::endl;
1878 }
1879 }
1880 delete it;
1881 TSynapse *synapse = nullptr;
1883 while ((synapse = (TSynapse *) it->Next())) {
1884 headerfile << " double synapse" << synapse << "();" << std::endl;
1885 sourcefile << "double " << classname << "::synapse"
1886 << synapse << "() {" << std::endl;
1887 sourcefile << " return (neuron" << synapse->GetPre()
1888 << "()*" << synapse->GetWeight() << ");" << std::endl;
1889 sourcefile << "}" << std::endl << std::endl;
1890 }
1891 delete it;
1892 headerfile << "};" << std::endl << std::endl;
1893 headerfile << "#endif // " << basefilename << "_h" << std::endl << std::endl;
1894 headerfile.close();
1895 sourcefile.close();
1896 std::cout << header << " and " << source << " created." << std::endl;
1897 }
1898 else if(lg == "FORTRAN") {
1899 TString implicit = " implicit double precision (a-h,n-z)\n";
1900 std::ofstream sigmoid("sigmoid.f");
1901 sigmoid << " double precision FUNCTION SIGMOID(X)" << std::endl
1902 << implicit
1903 << " IF(X.GT.37.) THEN" << std::endl
1904 << " SIGMOID = 1." << std::endl
1905 << " ELSE IF(X.LT.-709.) THEN" << std::endl
1906 << " SIGMOID = 0." << std::endl
1907 << " ELSE" << std::endl
1908 << " SIGMOID = 1./(1.+EXP(-X))" << std::endl
1909 << " ENDIF" << std::endl
1910 << " END" << std::endl;
1911 sigmoid.close();
1912 TString source = filename;
1913 source += ".f";
1914 std::ofstream sourcefile(source);
1915
1916 // Header
1917 sourcefile << " double precision function " << filename
1918 << "(x, index)" << std::endl;
1919 sourcefile << implicit;
1920 sourcefile << " double precision x(" <<
1921 fFirstLayer.GetEntriesFast() << ")" << std::endl << std::endl;
1922
1923 // Last layer
1924 sourcefile << "C --- Last Layer" << std::endl;
1925 TNeuron *neuron;
1927 Int_t idx = 0;
1928 TString ifelseif = " if (index.eq.";
1929 while ((neuron = (TNeuron *) it->Next())) {
1930 sourcefile << ifelseif.Data() << idx++ << ") then" << std::endl
1931 << " " << filename
1932 << "=neuron" << neuron << "(x);" << std::endl;
1933 ifelseif = " else if (index.eq.";
1934 }
1935 sourcefile << " else" << std::endl
1936 << " " << filename << "=0.d0" << std::endl
1937 << " endif" << std::endl;
1938 sourcefile << " end" << std::endl;
1939
1940 // Network
1941 sourcefile << "C --- First and Hidden layers" << std::endl;
1942 delete it;
1944 idx = 0;
1945 while ((neuron = (TNeuron *) it->Next())) {
1946 sourcefile << " double precision function neuron"
1947 << neuron << "(x)" << std::endl
1948 << implicit;
1949 sourcefile << " double precision x("
1950 << fFirstLayer.GetEntriesFast() << ")" << std::endl << std::endl;
1951 if (!neuron->GetPre(0)) {
1952 sourcefile << " neuron" << neuron
1953 << " = (x(" << idx+1 << ") - "
1954 << ((TNeuron *) fFirstLayer[idx])->GetNormalisation()[1]
1955 << "d0)/"
1956 << ((TNeuron *) fFirstLayer[idx])->GetNormalisation()[0]
1957 << "d0" << std::endl;
1958 idx++;
1959 } else {
1960 sourcefile << " neuron" << neuron
1961 << " = " << neuron->GetWeight() << "d0" << std::endl;
1962 TSynapse *syn;
1963 Int_t n = 0;
1964 while ((syn = neuron->GetPre(n++)))
1965 sourcefile << " neuron" << neuron
1966 << " = neuron" << neuron
1967 << " + synapse" << syn << "(x)" << std::endl;
1968 switch(neuron->GetType()) {
1969 case (TNeuron::kSigmoid):
1970 {
1971 sourcefile << " neuron" << neuron
1972 << "= (sigmoid(neuron" << neuron << ")*";
1973 break;
1974 }
1975 case (TNeuron::kLinear):
1976 {
1977 break;
1978 }
1979 case (TNeuron::kTanh):
1980 {
1981 sourcefile << " neuron" << neuron
1982 << "= (tanh(neuron" << neuron << ")*";
1983 break;
1984 }
1985 case (TNeuron::kGauss):
1986 {
1987 sourcefile << " neuron" << neuron
1988 << "= (exp(-neuron" << neuron << "*neuron"
1989 << neuron << "))*";
1990 break;
1991 }
1992 case (TNeuron::kSoftmax):
1993 {
1994 Int_t nn = 0;
1995 TNeuron* side = neuron->GetInLayer(nn++);
1996 sourcefile << " div = exp(neuron" << side << "())" << std::endl;
1997 while ((side = neuron->GetInLayer(nn++)))
1998 sourcefile << " div = div + exp(neuron" << side << "())" << std::endl;
1999 sourcefile << " neuron" << neuron ;
2000 sourcefile << "= (exp(neuron" << neuron << ") / div * ";
2001 break;
2002 }
2003 default:
2004 {
2005 sourcefile << " neuron " << neuron << "= 0.";
2006 }
2007 }
2008 sourcefile << neuron->GetNormalisation()[0] << "d0)+" ;
2009 sourcefile << neuron->GetNormalisation()[1] << "d0" << std::endl;
2010 }
2011 sourcefile << " end" << std::endl;
2012 }
2013 delete it;
2014
2015 // Synapses
2016 sourcefile << "C --- Synapses" << std::endl;
2017 TSynapse *synapse = nullptr;
2019 while ((synapse = (TSynapse *) it->Next())) {
2020 sourcefile << " double precision function " << "synapse"
2021 << synapse << "(x)\n" << implicit;
2022 sourcefile << " double precision x("
2023 << fFirstLayer.GetEntriesFast() << ")" << std::endl << std::endl;
2024 sourcefile << " synapse" << synapse
2025 << "=neuron" << synapse->GetPre()
2026 << "(x)*" << synapse->GetWeight() << "d0" << std::endl;
2027 sourcefile << " end" << std::endl << std::endl;
2028 }
2029 delete it;
2030 sourcefile.close();
2031 std::cout << source << " created." << std::endl;
2032 }
2033 else if(lg == "PYTHON") {
2034 TString classname = filename;
2035 TString pyfile = filename;
2036 pyfile += ".py";
2037 std::ofstream pythonfile(pyfile);
2038 pythonfile << "from math import exp" << std::endl << std::endl;
2039 pythonfile << "from math import tanh" << std::endl << std::endl;
2040 pythonfile << "class " << classname << ":" << std::endl;
2041 pythonfile << "\tdef value(self,index";
2042 for (i = 0; i < fFirstLayer.GetEntriesFast(); i++) {
2043 pythonfile << ",in" << i;
2044 }
2045 pythonfile << "):" << std::endl;
2046 for (i = 0; i < fFirstLayer.GetEntriesFast(); i++)
2047 pythonfile << "\t\tself.input" << i << " = (in" << i << " - "
2048 << ((TNeuron *) fFirstLayer[i])->GetNormalisation()[1] << ")/"
2049 << ((TNeuron *) fFirstLayer[i])->GetNormalisation()[0] << std::endl;
2050 TNeuron *neuron;
2052 Int_t idx = 0;
2053 while ((neuron = (TNeuron *) it->Next()))
2054 pythonfile << "\t\tif index==" << idx++
2055 << ": return self.neuron" << neuron << "();" << std::endl;
2056 pythonfile << "\t\treturn 0." << std::endl;
2057 delete it;
2059 idx = 0;
2060 while ((neuron = (TNeuron *) it->Next())) {
2061 pythonfile << "\tdef neuron" << neuron << "(self):" << std::endl;
2062 if (!neuron->GetPre(0))
2063 pythonfile << "\t\treturn self.input" << idx++ << std::endl;
2064 else {
2065 pythonfile << "\t\tinput = " << neuron->GetWeight() << std::endl;
2066 TSynapse *syn;
2067 Int_t n = 0;
2068 while ((syn = neuron->GetPre(n++)))
2069 pythonfile << "\t\tinput = input + self.synapse"
2070 << syn << "()" << std::endl;
2071 switch(neuron->GetType()) {
2072 case (TNeuron::kSigmoid):
2073 {
2074 pythonfile << "\t\tif input<-709. : return " << neuron->GetNormalisation()[1] << std::endl;
2075 pythonfile << "\t\treturn ((1/(1+exp(-input)))*";
2076 break;
2077 }
2078 case (TNeuron::kLinear):
2079 {
2080 pythonfile << "\t\treturn (input*";
2081 break;
2082 }
2083 case (TNeuron::kTanh):
2084 {
2085 pythonfile << "\t\treturn (tanh(input)*";
2086 break;
2087 }
2088 case (TNeuron::kGauss):
2089 {
2090 pythonfile << "\t\treturn (exp(-input*input)*";
2091 break;
2092 }
2093 case (TNeuron::kSoftmax):
2094 {
2095 pythonfile << "\t\treturn (exp(input) / (";
2096 Int_t nn = 0;
2097 TNeuron* side = neuron->GetInLayer(nn++);
2098 pythonfile << "exp(self.neuron" << side << "())";
2099 while ((side = neuron->GetInLayer(nn++)))
2100 pythonfile << " + exp(self.neuron" << side << "())";
2101 pythonfile << ") * ";
2102 break;
2103 }
2104 default:
2105 {
2106 pythonfile << "\t\treturn 0.";
2107 }
2108 }
2109 pythonfile << neuron->GetNormalisation()[0] << ")+" ;
2110 pythonfile << neuron->GetNormalisation()[1] << std::endl;
2111 }
2112 }
2113 delete it;
2114 TSynapse *synapse = nullptr;
2116 while ((synapse = (TSynapse *) it->Next())) {
2117 pythonfile << "\tdef synapse" << synapse << "(self):" << std::endl;
2118 pythonfile << "\t\treturn (self.neuron" << synapse->GetPre()
2119 << "()*" << synapse->GetWeight() << ")" << std::endl;
2120 }
2121 delete it;
2122 pythonfile.close();
2123 std::cout << pyfile << " created." << std::endl;
2124 }
2125}
2126
2127////////////////////////////////////////////////////////////////////////////////
2128/// Shuffle the Int_t index[n] in input.
2129///
2130/// Input:
2131/// - index: the array to shuffle
2132/// - n: the size of the array
2133///
2134/// Output:
2135/// - index: the shuffled indexes
2136///
2137/// This method is used for stochastic training
2138
2140{
2141 TTimeStamp ts;
2142 TRandom3 rnd(ts.GetSec());
2143 Int_t j, k;
2144 Int_t a = n - 1;
2145 for (Int_t i = 0; i < n; i++) {
2146 j = (Int_t) (rnd.Rndm() * a);
2147 k = index[j];
2148 index[j] = index[i];
2149 index[i] = k;
2150 }
2151 return;
2152}
2153
2154////////////////////////////////////////////////////////////////////////////////
2155/// One step for the stochastic method
2156/// buffer should contain the previous dw vector and will be updated
2157
2159{
2160 Int_t nEvents = fTraining->GetN();
2161 Int_t *index = new Int_t[nEvents];
2162 Int_t i,j,nentries;
2163 for (i = 0; i < nEvents; i++)
2164 index[i] = i;
2165 fEta *= fEtaDecay;
2166 Shuffle(index, nEvents);
2167 TNeuron *neuron;
2168 TSynapse *synapse;
2169 for (i = 0; i < nEvents; i++) {
2171 // First compute DeDw for all neurons: force calculation before
2172 // modifying the weights.
2174 for (j=0;j<nentries;j++) {
2175 neuron = (TNeuron *) fFirstLayer.UncheckedAt(j);
2176 neuron->GetDeDw();
2177 }
2178 Int_t cnt = 0;
2179 // Step for all neurons
2181 for (j=0;j<nentries;j++) {
2182 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2183 buffer[cnt] = (-fEta) * (neuron->GetDeDw() + fDelta)
2184 + fEpsilon * buffer[cnt];
2185 neuron->SetWeight(neuron->GetWeight() + buffer[cnt++]);
2186 }
2187 // Step for all synapses
2189 for (j=0;j<nentries;j++) {
2190 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2191 buffer[cnt] = (-fEta) * (synapse->GetDeDw() + fDelta)
2192 + fEpsilon * buffer[cnt];
2193 synapse->SetWeight(synapse->GetWeight() + buffer[cnt++]);
2194 }
2195 }
2196 delete[]index;
2197}
2198
2199////////////////////////////////////////////////////////////////////////////////
2200/// One step for the batch (stochastic) method.
2201/// DEDw should have been updated before calling this.
2202
2204{
2205 fEta *= fEtaDecay;
2206 Int_t cnt = 0;
2208 TNeuron *neuron = nullptr;
2209 // Step for all neurons
2210 while ((neuron = (TNeuron *) it->Next())) {
2211 buffer[cnt] = (-fEta) * (neuron->GetDEDw() + fDelta)
2212 + fEpsilon * buffer[cnt];
2213 neuron->SetWeight(neuron->GetWeight() + buffer[cnt++]);
2214 }
2215 delete it;
2217 TSynapse *synapse = nullptr;
2218 // Step for all synapses
2219 while ((synapse = (TSynapse *) it->Next())) {
2220 buffer[cnt] = (-fEta) * (synapse->GetDEDw() + fDelta)
2221 + fEpsilon * buffer[cnt];
2222 synapse->SetWeight(synapse->GetWeight() + buffer[cnt++]);
2223 }
2224 delete it;
2225}
2226
2227////////////////////////////////////////////////////////////////////////////////
2228/// Sets the weights to a point along a line
2229/// Weights are set to [origin + (dist * dir)].
2230
2232{
2233 Int_t idx = 0;
2234 TNeuron *neuron = nullptr;
2235 TSynapse *synapse = nullptr;
2237 while ((neuron = (TNeuron *) it->Next())) {
2238 neuron->SetWeight(origin[idx] + (dir[idx] * dist));
2239 idx++;
2240 }
2241 delete it;
2243 while ((synapse = (TSynapse *) it->Next())) {
2244 synapse->SetWeight(origin[idx] + (dir[idx] * dist));
2245 idx++;
2246 }
2247 delete it;
2248}
2249
2250////////////////////////////////////////////////////////////////////////////////
2251/// Sets the search direction to steepest descent.
2252
2254{
2255 Int_t idx = 0;
2256 TNeuron *neuron = nullptr;
2257 TSynapse *synapse = nullptr;
2259 while ((neuron = (TNeuron *) it->Next()))
2260 dir[idx++] = -neuron->GetDEDw();
2261 delete it;
2263 while ((synapse = (TSynapse *) it->Next()))
2264 dir[idx++] = -synapse->GetDEDw();
2265 delete it;
2266}
2267
2268////////////////////////////////////////////////////////////////////////////////
2269/// Search along the line defined by direction.
2270/// buffer is not used but is updated with the new dw
2271/// so that it can be used by a later stochastic step.
2272/// It returns true if the line search fails.
2273
2275{
2276 Int_t idx = 0;
2277 Int_t j,nentries;
2278 TNeuron *neuron = nullptr;
2279 TSynapse *synapse = nullptr;
2280 // store weights before line search
2281 Double_t *origin = new Double_t[fNetwork.GetEntriesFast() +
2284 for (j=0;j<nentries;j++) {
2285 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2286 origin[idx++] = neuron->GetWeight();
2287 }
2289 for (j=0;j<nentries;j++) {
2290 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2291 origin[idx++] = synapse->GetWeight();
2292 }
2293 // try to find a triplet (alpha1, alpha2, alpha3) such that
2294 // Error(alpha1)>Error(alpha2)<Error(alpha3)
2295 Double_t err1 = GetError(kTraining);
2296 Double_t alpha1 = 0.;
2297 Double_t alpha2 = fLastAlpha;
2298 if (alpha2 < 0.01)
2299 alpha2 = 0.01;
2300 if (alpha2 > 2.0)
2301 alpha2 = 2.0;
2302 Double_t alpha3 = alpha2;
2303 MLP_Line(origin, direction, alpha2);
2304 Double_t err2 = GetError(kTraining);
2305 Double_t err3 = err2;
2306 Bool_t bingo = false;
2307 Int_t icount;
2308 if (err1 > err2) {
2309 for (icount = 0; icount < 100; icount++) {
2310 alpha3 *= fTau;
2311 MLP_Line(origin, direction, alpha3);
2312 err3 = GetError(kTraining);
2313 if (err3 > err2) {
2314 bingo = true;
2315 break;
2316 }
2317 alpha1 = alpha2;
2318 err1 = err2;
2319 alpha2 = alpha3;
2320 err2 = err3;
2321 }
2322 if (!bingo) {
2323 MLP_Line(origin, direction, 0.);
2324 delete[]origin;
2325 return true;
2326 }
2327 } else {
2328 for (icount = 0; icount < 100; icount++) {
2329 alpha2 /= fTau;
2330 MLP_Line(origin, direction, alpha2);
2331 err2 = GetError(kTraining);
2332 if (err1 > err2) {
2333 bingo = true;
2334 break;
2335 }
2336 alpha3 = alpha2;
2337 err3 = err2;
2338 }
2339 if (!bingo) {
2340 MLP_Line(origin, direction, 0.);
2341 delete[]origin;
2342 fLastAlpha = 0.05;
2343 return true;
2344 }
2345 }
2346 // Sets the weights to the bottom of parabola
2347 fLastAlpha = 0.5 * (alpha1 + alpha3 -
2348 (err3 - err1) / ((err3 - err2) / (alpha3 - alpha2)
2349 - (err2 - err1) / (alpha2 - alpha1)));
2350 fLastAlpha = fLastAlpha < 10000 ? fLastAlpha : 10000;
2351 MLP_Line(origin, direction, fLastAlpha);
2353 // Stores weight changes (can be used by a later stochastic step)
2354 idx = 0;
2356 for (j=0;j<nentries;j++) {
2357 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2358 buffer[idx] = neuron->GetWeight() - origin[idx];
2359 idx++;
2360 }
2362 for (j=0;j<nentries;j++) {
2363 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2364 buffer[idx] = synapse->GetWeight() - origin[idx];
2365 idx++;
2366 }
2367 delete[]origin;
2368 return false;
2369}
2370
2371////////////////////////////////////////////////////////////////////////////////
2372/// Sets the search direction to conjugate gradient direction
2373/// beta should be:
2374///
2375/// \f$||g_{(t+1)}||^2 / ||g_{(t)}||^2\f$ (Fletcher-Reeves)
2376///
2377/// \f$g_{(t+1)} (g_{(t+1)}-g_{(t)}) / ||g_{(t)}||^2\f$ (Ribiere-Polak)
2378
2380{
2381 Int_t idx = 0;
2382 Int_t j,nentries;
2383 TNeuron *neuron = nullptr;
2384 TSynapse *synapse = nullptr;
2386 for (j=0;j<nentries;j++) {
2387 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2388 dir[idx] = -neuron->GetDEDw() + beta * dir[idx];
2389 idx++;
2390 }
2392 for (j=0;j<nentries;j++) {
2393 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2394 dir[idx] = -synapse->GetDEDw() + beta * dir[idx];
2395 idx++;
2396 }
2397}
2398
2399////////////////////////////////////////////////////////////////////////////////
2400/// Computes the hessian matrix using the BFGS update algorithm.
2401/// from gamma (g_{(t+1)}-g_{(t)}) and delta (w_{(t+1)}-w_{(t)}).
2402/// It returns true if such a direction could not be found
2403/// (if gamma and delta are orthogonal).
2404
2406{
2407 TMatrixD gd(gamma, TMatrixD::kTransposeMult, delta);
2408 if ((Double_t) gd[0][0] == 0.)
2409 return true;
2410 TMatrixD aHg(bfgsh, TMatrixD::kMult, gamma);
2411 TMatrixD tmp(gamma, TMatrixD::kTransposeMult, bfgsh);
2412 TMatrixD gHg(gamma, TMatrixD::kTransposeMult, aHg);
2413 Double_t a = 1 / (Double_t) gd[0][0];
2414 Double_t f = 1 + ((Double_t) gHg[0][0] * a);
2415 TMatrixD res( TMatrixD(delta, TMatrixD::kMult,
2417 res *= f;
2418 res -= (TMatrixD(delta, TMatrixD::kMult, tmp) +
2421 res *= a;
2422 bfgsh += res;
2423 return false;
2424}
2425
2426////////////////////////////////////////////////////////////////////////////////
2427/// Sets the gamma \f$(g_{(t+1)}-g_{(t)})\f$ and delta \f$(w_{(t+1)}-w_{(t)})\f$ vectors
2428/// Gamma is computed here, so ComputeDEDw cannot have been called before,
2429/// and delta is a direct translation of buffer into a TMatrixD.
2430
2432 Double_t * buffer)
2433{
2435 Int_t idx = 0;
2436 Int_t j,nentries;
2437 TNeuron *neuron = nullptr;
2438 TSynapse *synapse = nullptr;
2440 for (j=0;j<nentries;j++) {
2441 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2442 gamma[idx++][0] = -neuron->GetDEDw();
2443 }
2445 for (j=0;j<nentries;j++) {
2446 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2447 gamma[idx++][0] = -synapse->GetDEDw();
2448 }
2449 for (Int_t i = 0; i < els; i++)
2450 delta[i].Assign(buffer[i]);
2451 //delta.SetElements(buffer,"F");
2452 ComputeDEDw();
2453 idx = 0;
2455 for (j=0;j<nentries;j++) {
2456 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2457 gamma[idx++][0] += neuron->GetDEDw();
2458 }
2460 for (j=0;j<nentries;j++) {
2461 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2462 gamma[idx++][0] += synapse->GetDEDw();
2463 }
2464}
2465
2466////////////////////////////////////////////////////////////////////////////////
2467/// scalar product between gradient and direction
2468/// = derivative along direction
2469
2471{
2472 Int_t idx = 0;
2473 Int_t j,nentries;
2474 Double_t output = 0;
2475 TNeuron *neuron = nullptr;
2476 TSynapse *synapse = nullptr;
2478 for (j=0;j<nentries;j++) {
2479 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2480 output += neuron->GetDEDw() * dir[idx++];
2481 }
2483 for (j=0;j<nentries;j++) {
2484 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2485 output += synapse->GetDEDw() * dir[idx++];
2486 }
2487 return output;
2488}
2489
2490////////////////////////////////////////////////////////////////////////////////
2491/// Computes the direction for the BFGS algorithm as the product
2492/// between the Hessian estimate (bfgsh) and the dir.
2493
2495{
2497 TMatrixD dedw(els, 1);
2498 Int_t idx = 0;
2499 Int_t j,nentries;
2500 TNeuron *neuron = nullptr;
2501 TSynapse *synapse = nullptr;
2503 for (j=0;j<nentries;j++) {
2504 neuron = (TNeuron *) fNetwork.UncheckedAt(j);
2505 dedw[idx++][0] = neuron->GetDEDw();
2506 }
2508 for (j=0;j<nentries;j++) {
2509 synapse = (TSynapse *) fSynapses.UncheckedAt(j);
2510 dedw[idx++][0] = synapse->GetDEDw();
2511 }
2512 TMatrixD direction(bfgsh, TMatrixD::kMult, dedw);
2513 for (Int_t i = 0; i < els; i++)
2514 dir[i] = -direction[i][0];
2515 //direction.GetElements(dir,"F");
2516}
2517
2518////////////////////////////////////////////////////////////////////////////////
2519/// Draws the network structure.
2520/// Neurons are depicted by a blue disk, and synapses by
2521/// lines connecting neurons.
2522/// The line width is proportional to the weight.
2523
2525{
2526#define NeuronSize 2.5
2527
2528 Int_t nLayers = fStructure.CountChar(':')+1;
2529 Float_t xStep = 1./(nLayers+1.);
2530 Int_t layer;
2531 for(layer=0; layer< nLayers-1; layer++) {
2532 Float_t nNeurons_this = 0;
2533 if(layer==0) {
2535 nNeurons_this = input.CountChar(',')+1;
2536 }
2537 else {
2538 Int_t cnt=0;
2539 TString hidden = TString(fStructure(fStructure.First(':') + 1,fStructure.Last(':') - fStructure.First(':') - 1));
2540 Int_t beg = 0;
2541 Int_t end = hidden.Index(":", beg + 1);
2542 while (end != -1) {
2543 Int_t num = atoi(TString(hidden(beg, end - beg)).Data());
2544 cnt++;
2545 beg = end + 1;
2546 end = hidden.Index(":", beg + 1);
2547 if(layer==cnt) nNeurons_this = num;
2548 }
2549 Int_t num = atoi(TString(hidden(beg, hidden.Length() - beg)).Data());
2550 cnt++;
2551 if(layer==cnt) nNeurons_this = num;
2552 }
2553 Float_t nNeurons_next = 0;
2554 if(layer==nLayers-2) {
2556 nNeurons_next = output.CountChar(',')+1;
2557 }
2558 else {
2559 Int_t cnt=0;
2560 TString hidden = TString(fStructure(fStructure.First(':') + 1,fStructure.Last(':') - fStructure.First(':') - 1));
2561 Int_t beg = 0;
2562 Int_t end = hidden.Index(":", beg + 1);
2563 while (end != -1) {
2564 Int_t num = atoi(TString(hidden(beg, end - beg)).Data());
2565 cnt++;
2566 beg = end + 1;
2567 end = hidden.Index(":", beg + 1);
2568 if(layer+1==cnt) nNeurons_next = num;
2569 }
2570 Int_t num = atoi(TString(hidden(beg, hidden.Length() - beg)).Data());
2571 cnt++;
2572 if(layer+1==cnt) nNeurons_next = num;
2573 }
2574 Float_t yStep_this = 1./(nNeurons_this+1.);
2575 Float_t yStep_next = 1./(nNeurons_next+1.);
2577 TSynapse *theSynapse = nullptr;
2578 Float_t maxWeight = 0;
2579 while ((theSynapse = (TSynapse *) it->Next()))
2580 maxWeight = maxWeight < theSynapse->GetWeight() ? theSynapse->GetWeight() : maxWeight;
2581 delete it;
2583 for(Int_t neuron1=0; neuron1<nNeurons_this; neuron1++) {
2584 for(Int_t neuron2=0; neuron2<nNeurons_next; neuron2++) {
2585 TLine* synapse = new TLine(xStep*(layer+1),yStep_this*(neuron1+1),xStep*(layer+2),yStep_next*(neuron2+1));
2586 synapse->Draw();
2587 theSynapse = (TSynapse *) it->Next();
2588 if (!theSynapse) continue;
2589 synapse->SetLineWidth(Int_t((theSynapse->GetWeight()/maxWeight)*10.));
2590 synapse->SetLineStyle(1);
2591 if(((TMath::Abs(theSynapse->GetWeight())/maxWeight)*10.)<0.5) synapse->SetLineStyle(2);
2592 if(((TMath::Abs(theSynapse->GetWeight())/maxWeight)*10.)<0.25) synapse->SetLineStyle(3);
2593 }
2594 }
2595 delete it;
2596 }
2597 for(layer=0; layer< nLayers; layer++) {
2598 Float_t nNeurons = 0;
2599 if(layer==0) {
2601 nNeurons = input.CountChar(',')+1;
2602 }
2603 else if(layer==nLayers-1) {
2605 nNeurons = output.CountChar(',')+1;
2606 }
2607 else {
2608 Int_t cnt=0;
2609 TString hidden = TString(fStructure(fStructure.First(':') + 1,fStructure.Last(':') - fStructure.First(':') - 1));
2610 Int_t beg = 0;
2611 Int_t end = hidden.Index(":", beg + 1);
2612 while (end != -1) {
2613 Int_t num = atoi(TString(hidden(beg, end - beg)).Data());
2614 cnt++;
2615 beg = end + 1;
2616 end = hidden.Index(":", beg + 1);
2617 if(layer==cnt) nNeurons = num;
2618 }
2619 Int_t num = atoi(TString(hidden(beg, hidden.Length() - beg)).Data());
2620 cnt++;
2621 if(layer==cnt) nNeurons = num;
2622 }
2623 Float_t yStep = 1./(nNeurons+1.);
2624 for(Int_t neuron=0; neuron<nNeurons; neuron++) {
2625 TMarker* m = new TMarker(xStep*(layer+1),yStep*(neuron+1),20);
2626 m->SetMarkerColor(4);
2628 m->Draw();
2629 }
2630 }
2631 const TString input = TString(fStructure(0, fStructure.First(':')));
2632 const TObjArray *inpL = input.Tokenize(" ,");
2633 const Int_t nrItems = inpL->GetLast()+1;
2634 Float_t yStep = 1./(nrItems+1);
2635 for (Int_t item = 0; item < nrItems; item++) {
2636 const TString brName = ((TObjString *)inpL->At(item))->GetString();
2637 TText* label = new TText(0.5*xStep,yStep*(item+1),brName.Data());
2638 label->Draw();
2639 }
2640 delete inpL;
2641
2642 Int_t numOutNodes=fLastLayer.GetEntriesFast();
2643 yStep=1./(numOutNodes+1);
2644 for (Int_t outnode=0; outnode<numOutNodes; outnode++) {
2645 TNeuron* neuron=(TNeuron*)fLastLayer[outnode];
2646 if (neuron && neuron->GetName()) {
2647 TText* label = new TText(xStep*nLayers,
2648 yStep*(outnode+1),
2649 neuron->GetName());
2650 label->Draw();
2651 }
2652 }
2653}
#define f(i)
Definition RSha256.hxx:104
#define a(i)
Definition RSha256.hxx:99
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Ssiz_t kNPOS
Definition RtypesCore.h:117
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
#define gDirectory
Definition TDirectory.h:384
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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 filename
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 target
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
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
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void reg
char name[80]
Definition TGX11.cxx:110
int nentries
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
#define NeuronSize
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
#define gPad
static void SetMaxima(Int_t maxop=1000, Int_t maxpar=1000, Int_t maxconst=1000)
static function to set the maximum value of 3 parameters
static void GetMaxima(Int_t &maxop, Int_t &maxpar, Int_t &maxconst)
static function to get the maximum value of 3 parameters -maxop : maximum number of operations -maxpa...
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
virtual void SetLeftMargin(Float_t leftmargin)
Set Pad left margin in fraction of the pad width.
Definition TAttPad.cxx:109
void SetDecimals(Bool_t dot=kTRUE)
Sets the decimals flag By default, blank characters are stripped, and then the label is correctly ali...
Definition TAxis.h:213
virtual void UnZoom()
Reset first & last bin to the full range.
Definition TAxis.cxx:1269
The Canvas class.
Definition TCanvas.h:23
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition TClass.cxx:3037
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
<div class="legacybox"><h2>Legacy Code</h2> TEventList is a legacy interface: there will be no bug fi...
Definition TEventList.h:31
virtual Long64_t GetEntry(Int_t index) const
Return value of entry at index in the list.
virtual Int_t GetN() const
Definition TEventList.h:56
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2342
Double_t * GetY() const
Definition TGraph.h:140
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:670
void Reset(Option_t *option="") override
Reset.
Definition TH1.cxx:10512
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
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:357
void Reset(Option_t *option="") override
Reset this histogram: contents, errors, etc.
Definition TH2.cxx:4246
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:393
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
Use the TLine constructor to create a simple line.
Definition TLine.h:22
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
Manages Markers.
Definition TMarker.h:22
void Draw(Option_t *option="") override
Draw this marker with its current attributes.
Definition TMarker.cxx:199
virtual TMatrixTBase< Element > & UnitMatrix()
Make a unit matrix (matrix need not be a square one).
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
TList * GetListOfGraphs() const
Definition TMultiGraph.h:68
virtual void Add(TGraph *graph, Option_t *chopt="")
Add a new graph to the list of graphs.
void Draw(Option_t *chopt="") override
Draw this multigraph with its current attributes.
TAxis * GetYaxis()
Get y axis of the graph.
TAxis * GetXaxis()
Get x axis of the graph.
This class describes a neural network.
TTreeFormula * fEventWeight
! formula representing the event weight
void BuildOneHiddenLayer(const TString &sNumNodes, Int_t &layer, Int_t &prevStart, Int_t &prevStop, Bool_t lastLayer)
Builds a hidden layer, updates the number of layers.
void SteepestDir(Double_t *)
Sets the search direction to steepest descent.
void BuildNetwork()
Instantiates the network from the description.
TObjArray fNetwork
Collection of all the neurons in the 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.
TEventList * fTest
! EventList defining the events in the test dataset
bool GetBFGSH(TMatrixD &, TMatrixD &, TMatrixD &)
Computes the hessian matrix using the BFGS update algorithm.
void BuildHiddenLayers(TString &)
Builds hidden layers.
void BuildFirstLayer(TString &)
Instantiates the neurons in input Inputs are normalised and the type is set to kOff (simple forward o...
void SetTau(Double_t tau)
Sets Tau - used in line search (look at the constructor for the complete description of learning meth...
TMultiLayerPerceptron()
Default constructor.
Double_t GetSumSquareError() const
Error on the output for a given event.
void ConjugateGradientsDir(Double_t *, Double_t)
Sets the search direction to conjugate gradient direction beta should be:
Double_t fTau
! Tau - used in line search - Default=3.
TTree * fData
! pointer to the tree used as datasource
Double_t Result(Int_t event, Int_t index=0) const
Computes the output for a given event.
void SetGammaDelta(TMatrixD &, TMatrixD &, Double_t *)
Sets the gamma and delta vectors Gamma is computed here, so ComputeDEDw cannot have been called bef...
TEventList * fTraining
! EventList defining the events in the training dataset
TString fStructure
String containing the network structure.
Int_t fReset
! number of epochs between two resets of the search direction to the steepest descent - Default=50
Bool_t LoadWeights(Option_t *filename="")
Loads the weights from a text file conforming to the format defined by DumpWeights.
void MLP_Batch(Double_t *)
One step for the batch (stochastic) method.
TNeuron::ENeuronType fOutType
Type of output neurons.
Double_t fCurrentTreeWeight
! weight of the current tree in a chain
ELearningMethod fLearningMethod
! The Learning Method
Double_t fLastAlpha
! internal parameter used in line search
Int_t fCurrentTree
! index of the current tree in a chain
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++ ,...
Double_t fEpsilon
! Epsilon - used in stochastic minimisation - Default=0.
void Train(Int_t nEpoch, Option_t *option="text", Double_t minE=0)
Train the network.
TNeuron::ENeuronType GetType() const
void BFGSDir(TMatrixD &, Double_t *)
Computes the direction for the BFGS algorithm as the product between the Hessian estimate (bfgsh) and...
void SetTestDataSet(TEventList *test)
Sets the Test dataset.
Bool_t fTrainingOwner
! internal flag whether one has to delete fTraining or not
void SetLearningMethod(TMultiLayerPerceptron::ELearningMethod method)
Sets the learning method.
void SetTrainingDataSet(TEventList *train)
Sets the Training dataset.
void BuildLastLayer(TString &, Int_t)
Builds the output layer Neurons are linear combinations of input, by default.
Double_t fDelta
! Delta - used in stochastic minimisation - Default=0.
TTreeFormulaManager * fManager
! TTreeFormulaManager for the weight and neurons
void Randomize() const
Randomize the weights.
Bool_t LineSearch(Double_t *, Double_t *)
Search along the line defined by direction.
void ExpandStructure()
Expand the structure of the first layer.
Double_t fEta
! Eta - used in stochastic minimisation - Default=0.1
Double_t GetError(Int_t event) const
Error on the output for a given event.
Double_t fEtaDecay
! EtaDecay - Eta *= EtaDecay at each epoch - Default=1.
void SetEtaDecay(Double_t ed)
Sets EtaDecay - Eta *= EtaDecay at each epoch (look at the constructor for the complete description o...
void AttachData()
Connects the TTree to Neurons in input and output layers.
void SetData(TTree *)
Set the data source.
void SetEventWeight(const char *)
Set the event weight.
Bool_t DumpWeights(Option_t *filename="-") const
Dumps the weights to a text file.
TString fWeight
String containing the event weight.
void SetDelta(Double_t delta)
Sets Delta - used in stochastic minimisation (look at the constructor for the complete description of...
~TMultiLayerPerceptron() override
Destructor.
Double_t GetCrossEntropy() const
Cross entropy error for a softmax output neuron, for a given event.
void SetReset(Int_t reset)
Sets number of epochs between two resets of the search direction to the steepest descent.
Bool_t fTestOwner
! internal flag whether one has to delete fTest or not
void Shuffle(Int_t *, Int_t) const
Shuffle the Int_t index[n] in input.
Double_t DerivDir(Double_t *)
scalar product between gradient and direction = derivative along direction
void MLP_Stochastic(Double_t *)
One step for the stochastic method buffer should contain the previous dw vector and will be updated.
void Draw(Option_t *option="") override
Draws the network structure.
TObjArray fSynapses
Collection of all the synapses in the network.
void MLP_Line(Double_t *, Double_t *, Double_t)
Sets the weights to a point along a line Weights are set to [origin + (dist * dir)].
TNeuron::ENeuronType fType
Type of hidden neurons.
TObjArray fLastLayer
Collection of the output neurons; subset of fNetwork.
TString fextD
String containing the derivative name.
void ComputeDEDw() const
Compute the DEDw = sum on all training events of dedw for each weight normalized by the number of eve...
Double_t GetCrossEntropyBinary() const
Cross entropy error for sigmoid output neurons, for a given event.
void DrawResult(Int_t index=0, Option_t *option="test") const
Draws the neural net output It produces an histogram with the output for the two datasets.
void SetEta(Double_t eta)
Sets Eta - used in stochastic minimisation (look at the constructor for the complete description of l...
TObjArray fFirstLayer
Collection of the input neurons; subset of fNetwork.
void GetEntry(Int_t) const
Load an entry into the network.
void SetEpsilon(Double_t eps)
Sets Epsilon - used in stochastic minimisation (look at the constructor for the complete description ...
TString fextF
String containing the function name.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
This class describes an elementary neuron, which is the basic element for a Neural Network.
Definition TNeuron.h:25
Double_t GetWeight() const
Definition TNeuron.h:48
void SetWeight(Double_t w)
Sets the neuron weight to w.
Definition TNeuron.cxx:1144
Double_t GetDEDw() const
Definition TNeuron.h:53
Double_t GetValue() const
Computes the output using the appropriate function and all the weighted inputs, or uses the branch as...
Definition TNeuron.cxx:944
void SetDEDw(Double_t in)
Sets the derivative of the total error wrt the neuron weight.
Definition TNeuron.cxx:1164
Double_t GetDeDw() const
Computes the derivative of the error wrt the neuron weight.
Definition TNeuron.cxx:1080
TNeuron * GetInLayer(Int_t n) const
Definition TNeuron.h:37
Double_t GetError() const
Computes the error for output neurons.
Definition TNeuron.cxx:1059
TTreeFormula * UseBranch(TTree *, const char *)
Sets a formula that can be used to make the neuron an input.
Definition TNeuron.cxx:874
TSynapse * GetPre(Int_t n) const
Definition TNeuron.h:35
void ForceExternalValue(Double_t value)
Uses the branch type to force an external value.
Definition TNeuron.cxx:1121
Double_t GetTarget() const
Computes the normalized target pattern for output neurons.
Definition TNeuron.cxx:1070
const Double_t * GetNormalisation() const
Definition TNeuron.h:50
ENeuronType GetType() const
Returns the neuron type.
Definition TNeuron.cxx:863
void SetNewEvent() const
Inform the neuron that inputs of the network have changed, so that the buffered values have to be rec...
Definition TNeuron.cxx:1153
void SetNormalisation(Double_t mean, Double_t RMS)
Sets the normalization variables.
Definition TNeuron.cxx:1133
void AddInLayer(TNeuron *)
Tells a neuron which neurons form its layer (including itself).
Definition TNeuron.cxx:853
ENeuronType
Definition TNeuron.h:29
@ kLinear
Definition TNeuron.h:29
@ kSigmoid
Definition TNeuron.h:29
@ kExternal
Definition TNeuron.h:29
@ kSoftmax
Definition TNeuron.h:29
@ kGauss
Definition TNeuron.h:29
@ kOff
Definition TNeuron.h:29
@ kTanh
Definition TNeuron.h:29
Iterator of object array.
Definition TObjArray.h:117
TObject * Next() override
Return next object in array. Returns 0 when no more objects in array.
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
TIterator * MakeIterator(Bool_t dir=kIterForward) const override
Returns an array iterator.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
void AddLast(TObject *obj) override
Add object in the next empty slot in the array.
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
Int_t GetLast() const override
Return index of last object in array.
Collectable string class.
Definition TObjString.h:28
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:979
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:280
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:967
Random number generator class based on M.
Definition TRandom3.h:27
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom3.cxx:99
Regular expression class.
Definition TRegexp.h:31
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:417
void ToLower()
Change string to lower-case.
Definition TString.cxx:1182
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Return true if string ends with the specified string.
Definition TString.cxx:2244
Ssiz_t First(char c) const
Find first occurrence of a character c.
Definition TString.cxx:538
const char * Data() const
Definition TString.h:376
Bool_t IsAlpha() const
Returns true if all characters in string are alphabetic.
Definition TString.cxx:1798
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition TString.cxx:931
void ToUpper()
Change string to upper case.
Definition TString.cxx:1195
Int_t CountChar(Int_t c) const
Return number of times character c occurs in the string.
Definition TString.cxx:515
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
Bool_t IsAlnum() const
Returns true if all characters in string are alphanumeric.
Definition TString.cxx:1813
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
This is a simple weighted bidirectional connection between two neurons.
Definition TSynapse.h:20
Double_t GetDeDw() const
Computes the derivative of the error wrt the synapse weight.
Definition TSynapse.cxx:89
void SetWeight(Double_t w)
Sets the weight of the synapse.
Definition TSynapse.cxx:102
Double_t GetWeight() const
Definition TSynapse.h:30
void SetDEDw(Double_t in)
Sets the derivative of the total error wrt the synapse weight.
Definition TSynapse.cxx:110
TNeuron * GetPre() const
Definition TSynapse.h:27
Double_t GetDEDw() const
Definition TSynapse.h:34
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition TSystem.cxx:1857
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:416
Base class for several text objects.
Definition TText.h:22
The TTimeStamp encapsulates seconds and ns since EPOCH.
Definition TTimeStamp.h:45
time_t GetSec() const
Definition TTimeStamp.h:109
Used to coordinate one or more TTreeFormula objects.
bool Notify() override
This method must be overridden to handle object notification (the base implementation is no-op).
virtual void Add(TTreeFormula *)
Add a new formula to the list of formulas managed The manager of the formula will be changed and the ...
virtual void Remove(TTreeFormula *)
Remove a formula from this manager.
Used to pass a selection expression to the Tree drawing routine.
T EvalInstance(Int_t i=0, const char *stringStack[]=nullptr)
Evaluate this treeformula.
A TTree represents a columnar dataset.
Definition TTree.h:79
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 Double_t GetWeight() const
Definition TTree.h:584
void Draw(Option_t *opt) override
Default Draw method for all objects.
Definition TTree.h:431
virtual Long64_t GetEntries() const
Definition TTree.h:463
virtual Int_t GetTreeNumber() const
Definition TTree.h:559
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition TVirtualPad.h:51
virtual void Modified(Bool_t flag=1)=0
virtual void Update()=0
const Int_t n
Definition legend1.C:16
Bool_t IsNaN(Double_t x)
Definition TMath.h:892
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:756
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
TCanvas * slash()
Definition slash.C:1
TMarker m
Definition textangle.C:8
static void output()