Logo ROOT   6.14/05
Reference Guide
TMultiLayerPerceptron.h
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 #ifndef ROOT_TMultiLayerPerceptron
13 #define ROOT_TMultiLayerPerceptron
14 
15 #include "TObject.h"
16 #include "TString.h"
17 #include "TObjArray.h"
18 #include "TMatrixD.h"
19 #include "TNeuron.h"
20 
21 class TTree;
22 class TEventList;
23 class TTreeFormula;
25 
26 //____________________________________________________________________
27 //
28 // TMultiLayerPerceptron
29 //
30 // This class decribes a Neural network.
31 // There are facilities to train the network and use the output.
32 //
33 // The input layer is made of inactive neurons (returning the
34 // normalized input), hidden layers are made of sigmoids and output
35 // neurons are linear.
36 //
37 // The basic input is a TTree and two (training and test) TEventLists.
38 // For classification jobs, a branch (maybe in a TFriend) must contain
39 // the expected output.
40 // 6 learning methods are available: kStochastic, kBatch,
41 // kSteepestDescent, kRibierePolak, kFletcherReeves and kBFGS.
42 //
43 // This implementation is *inspired* from the mlpfit package from
44 // J.Schwindling et al.
45 //
46 //____________________________________________________________________
47 
49  friend class TMLPAnalyzer;
50 
51  public:
56  TMultiLayerPerceptron(const char* layout, TTree* data = 0,
57  const char* training = "Entry$%2==0",
58  const char* test = "",
60  const char* extF = "", const char* extD = "");
61  TMultiLayerPerceptron(const char* layout,
62  const char* weight, TTree* data = 0,
63  const char* training = "Entry$%2==0",
64  const char* test = "",
66  const char* extF = "", const char* extD = "");
67  TMultiLayerPerceptron(const char* layout, TTree* data,
68  TEventList* training,
71  const char* extF = "", const char* extD = "");
72  TMultiLayerPerceptron(const char* layout,
73  const char* weight, TTree* data,
74  TEventList* training,
75  TEventList* test,
77  const char* extF = "", const char* extD = "");
78  virtual ~TMultiLayerPerceptron();
79  void SetData(TTree*);
80  void SetTrainingDataSet(TEventList* train);
81  void SetTestDataSet(TEventList* test);
82  void SetTrainingDataSet(const char* train);
83  void SetTestDataSet(const char* test);
85  void SetEventWeight(const char*);
86  void Train(Int_t nEpoch, Option_t* option = "text", Double_t minE=0);
87  Double_t Result(Int_t event, Int_t index = 0) const;
88  Double_t GetError(Int_t event) const;
90  void ComputeDEDw() const;
91  void Randomize() const;
92  void SetEta(Double_t eta);
93  void SetEpsilon(Double_t eps);
94  void SetDelta(Double_t delta);
95  void SetEtaDecay(Double_t ed);
96  void SetTau(Double_t tau);
97  void SetReset(Int_t reset);
98  inline Double_t GetEta() const { return fEta; }
99  inline Double_t GetEpsilon() const { return fEpsilon; }
100  inline Double_t GetDelta() const { return fDelta; }
101  inline Double_t GetEtaDecay() const { return fEtaDecay; }
103  inline Double_t GetTau() const { return fTau; }
104  inline Int_t GetReset() const { return fReset; }
105  inline TString GetStructure() const { return fStructure; }
106  inline TNeuron::ENeuronType GetType() const { return fType; }
107  void DrawResult(Int_t index = 0, Option_t* option = "test") const;
108  Bool_t DumpWeights(Option_t* filename = "-") const;
109  Bool_t LoadWeights(Option_t* filename = "");
110  Double_t Evaluate(Int_t index, Double_t* params) const;
111  void Export(Option_t* filename = "NNfunction", Option_t* language = "C++") const;
112  virtual void Draw(Option_t *option="");
113 
114  protected:
115  void AttachData();
116  void BuildNetwork();
117  void GetEntry(Int_t) const;
118  // it's a choice not to force learning function being const, even if possible
119  void MLP_Stochastic(Double_t*);
120  void MLP_Batch(Double_t*);
122  void SteepestDir(Double_t*);
125  bool GetBFGSH(TMatrixD&, TMatrixD &, TMatrixD&);
126  void BFGSDir(TMatrixD&, Double_t*);
129  Double_t GetCrossEntropy() const;
130  Double_t GetSumSquareError() const;
131 
132  private:
133  TMultiLayerPerceptron(const TMultiLayerPerceptron&); // Not implemented
134  TMultiLayerPerceptron& operator=(const TMultiLayerPerceptron&); // Not implemented
135  void ExpandStructure();
136  void BuildFirstLayer(TString&);
137  void BuildHiddenLayers(TString&);
138  void BuildOneHiddenLayer(const TString& sNumNodes, Int_t& layer,
139  Int_t& prevStart, Int_t& prevStop,
140  Bool_t lastLayer);
141  void BuildLastLayer(TString&, Int_t);
142  void Shuffle(Int_t*, Int_t) const;
144 
145  TTree* fData; //! pointer to the tree used as datasource
146  Int_t fCurrentTree; //! index of the current tree in a chain
147  Double_t fCurrentTreeWeight; //! weight of the current tree in a chain
148  TObjArray fNetwork; // Collection of all the neurons in the network
149  TObjArray fFirstLayer; // Collection of the input neurons; subset of fNetwork
150  TObjArray fLastLayer; // Collection of the output neurons; subset of fNetwork
151  TObjArray fSynapses; // Collection of all the synapses in the network
152  TString fStructure; // String containing the network structure
153  TString fWeight; // String containing the event weight
154  TNeuron::ENeuronType fType; // Type of hidden neurons
155  TNeuron::ENeuronType fOutType; // Type of output neurons
156  TString fextF; // String containing the function name
157  TString fextD; // String containing the derivative name
158  TEventList *fTraining; //! EventList defining the events in the training dataset
159  TEventList *fTest; //! EventList defining the events in the test dataset
160  ELearningMethod fLearningMethod; //! The Learning Method
161  TTreeFormula* fEventWeight; //! formula representing the event weight
162  TTreeFormulaManager* fManager; //! TTreeFormulaManager for the weight and neurons
163  Double_t fEta; //! Eta - used in stochastic minimisation - Default=0.1
164  Double_t fEpsilon; //! Epsilon - used in stochastic minimisation - Default=0.
165  Double_t fDelta; //! Delta - used in stochastic minimisation - Default=0.
166  Double_t fEtaDecay; //! EtaDecay - Eta *= EtaDecay at each epoch - Default=1.
167  Double_t fTau; //! Tau - used in line search - Default=3.
168  Double_t fLastAlpha; //! internal parameter used in line search
169  Int_t fReset; //! number of epochs between two resets of the search direction to the steepest descent - Default=50
170  Bool_t fTrainingOwner; //! internal flag whether one has to delete fTraining or not
171  Bool_t fTestOwner; //! internal flag whether one has to delete fTest or not
172  ClassDef(TMultiLayerPerceptron, 4) // a Neural Network
173 };
174 
175 #endif
Double_t DerivDir(Double_t *)
scalar product between gradient and direction = derivative along direction
Bool_t LineSearch(Double_t *, Double_t *)
Search along the line defined by direction.
void SteepestDir(Double_t *)
Sets the search direction to steepest descent.
An array of TObjects.
Definition: TObjArray.h:37
void Randomize() const
Randomize the weights.
void BuildHiddenLayers(TString &)
Builds hidden layers.
void BFGSDir(TMatrixD &, Double_t *)
Computes the direction for the BFGS algorithm as the product between the Hessian estimate (bfgsh) and...
Double_t fEpsilon
Eta - used in stochastic minimisation - Default=0.1.
TTreeFormulaManager * fManager
formula representing the event weight
const char Option_t
Definition: RtypesCore.h:62
void SetTestDataSet(TEventList *test)
Sets the Test dataset.
Double_t fCurrentTreeWeight
index of the current tree in a chain
void SetEpsilon(Double_t eps)
Sets Epsilon - used in stochastic minimisation (look at the constructor for the complete description ...
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 Evaluate(Int_t index, Double_t *params) const
Returns the Neural Net for a given set of input parameters #parameters must equal #input neurons...
Definition: test.py:1
ENeuronType
Definition: TNeuron.h:48
void AttachData()
Connects the TTree to Neurons in input and output layers.
Basic string class.
Definition: TString.h:131
void SetEtaDecay(Double_t ed)
Sets EtaDecay - Eta *= EtaDecay at each epoch (look at the constructor for the complete description o...
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
void SetEta(Double_t eta)
Sets Eta - used in stochastic minimisation (look at the constructor for the complete description of l...
TMultiLayerPerceptron::ELearningMethod GetLearningMethod() const
void SetData(TTree *)
Set the data source.
void Shuffle(Int_t *, Int_t) const
Shuffle the Int_t index[n] in input.
Double_t GetSumSquareError() const
Error on the output for a given event.
TObjArray fNetwork
weight of the current tree in a chain
void ConjugateGradientsDir(Double_t *, Double_t)
Sets the search direction to conjugate gradient direction beta should be: ||g_{(t+1)}||^2 / ||g_{(t)}...
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.
Used to coordinate one or more TTreeFormula objects.
#define ClassDef(name, id)
Definition: Rtypes.h:320
bool GetBFGSH(TMatrixD &, TMatrixD &, TMatrixD &)
Computes the hessian matrix using the BFGS update algorithm.
Double_t fEta
TTreeFormulaManager for the weight and neurons.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
TMultiLayerPerceptron & operator=(const TMultiLayerPerceptron &)
virtual ~TMultiLayerPerceptron()
Destructor.
void SetReset(Int_t reset)
Sets number of epochs between two resets of the search direction to the steepest descent.
Bool_t DumpWeights(Option_t *filename="-") const
Dumps the weights to a text file.
Int_t fCurrentTree
pointer to the tree used as datasource
TNeuron::ENeuronType fType
Used to pass a selection expression to the Tree drawing routine.
Definition: TTreeFormula.h:58
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 Train(Int_t nEpoch, Option_t *option="text", Double_t minE=0)
Train the network.
TNeuron::ENeuronType fOutType
virtual void Draw(Option_t *option="")
Draws the network structure.
Bool_t fTestOwner
internal flag whether one has to delete fTraining or not
Double_t GetError(Int_t event) const
Error on the output for a given event.
A TEventList object is a list of selected events (entries) in a TTree.
Definition: TEventList.h:31
Bool_t LoadWeights(Option_t *filename="")
Loads the weights from a text file conforming to the format defined by DumpWeights.
Int_t fReset
internal parameter used in line search
void SetTrainingDataSet(TEventList *train)
Sets the Training dataset.
void SetGammaDelta(TMatrixD &, TMatrixD &, Double_t *)
Sets the gamma (g_{(t+1)}-g_{(t)}) and delta (w_{(t+1)}-w_{(t)}) vectors Gamma is computed here...
double Double_t
Definition: RtypesCore.h:55
void ComputeDEDw() const
Compute the DEDw = sum on all training events of dedw for each weight normalized by the number of eve...
Double_t fEtaDecay
Delta - used in stochastic minimisation - Default=0.
Double_t fDelta
Epsilon - used in stochastic minimisation - Default=0.
int type
Definition: TGX11.cxx:120
void SetDelta(Double_t delta)
Sets Delta - used in stochastic minimisation (look at the constructor for the complete description of...
Double_t fTau
EtaDecay - Eta *= EtaDecay at each epoch - Default=1.
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)].
TEventList * fTest
EventList defining the events in the training dataset.
void MLP_Batch(Double_t *)
One step for the batch (stochastic) method.
Mother of all ROOT objects.
Definition: TObject.h:37
ELearningMethod fLearningMethod
EventList defining the events in the test dataset.
void MLP_Stochastic(Double_t *)
One step for the stochastic method buffer should contain the previous dw vector and will be updated...
Bool_t fTrainingOwner
number of epochs between two resets of the search direction to the steepest descent - Default=50 ...
Double_t GetCrossEntropyBinary() const
Cross entropy error for sigmoid output neurons, for a given event.
void ExpandStructure()
Expand the structure of the first layer.
TMultiLayerPerceptron()
Default constructor.
A TTree object has a header with a name and a title.
Definition: TTree.h:70
void BuildFirstLayer(TString &)
Instanciates the neurons in input Inputs are normalised and the type is set to kOff (simple forward o...
void SetEventWeight(const char *)
Set the event weight.
void GetEntry(Int_t) const
Load an entry into the network.
Double_t fLastAlpha
Tau - used in line search - Default=3.
Double_t Result(Int_t event, Int_t index=0) const
Computes the output for a given event.
Double_t GetCrossEntropy() const
Cross entropy error for a softmax output neuron, for a given event.
void BuildLastLayer(TString &, Int_t)
Builds the output layer Neurons are linear combinations of input, by defaul.
void SetLearningMethod(TMultiLayerPerceptron::ELearningMethod method)
Sets the learning method.
void SetTau(Double_t tau)
Sets Tau - used in line search (look at the constructor for the complete description of learning meth...
void BuildNetwork()
Instanciates the network from the description.
TNeuron::ENeuronType GetType() const
TTreeFormula * fEventWeight
The Learning Method.