#ifndef ROOT_TMultiLayerPerceptron
#define ROOT_TMultiLayerPerceptron
#ifndef ROOT_TObject
#include "TObject.h"
#endif
#ifndef ROOT_TString
#include "TString.h"
#endif
#ifndef ROOT_TObjArray
#include "TObjArray.h"
#endif
#ifndef ROOT_TMatrixD
#include "TMatrixD.h"
#endif
#ifndef ROOT_TNeuron
#include "TNeuron.h"
#endif
class TTree;
class TEventList;
class TTreeFormula;
class TTreeFormulaManager;
class TMultiLayerPerceptron : public TObject {
friend class TMLPAnalyzer;
public:
enum ELearningMethod { kStochastic, kBatch, kSteepestDescent,
kRibierePolak, kFletcherReeves, kBFGS };
enum EDataSet { kTraining, kTest };
TMultiLayerPerceptron();
TMultiLayerPerceptron(const char* layout, TTree* data = 0,
const char* training = "Entry$%2==0",
const char* test = "",
TNeuron::ENeuronType type = TNeuron::kSigmoid,
const char* extF = "", const char* extD = "");
TMultiLayerPerceptron(const char* layout,
const char* weight, TTree* data = 0,
const char* training = "Entry$%2==0",
const char* test = "",
TNeuron::ENeuronType type = TNeuron::kSigmoid,
const char* extF = "", const char* extD = "");
TMultiLayerPerceptron(const char* layout, TTree* data,
TEventList* training,
TEventList* test,
TNeuron::ENeuronType type = TNeuron::kSigmoid,
const char* extF = "", const char* extD = "");
TMultiLayerPerceptron(const char* layout,
const char* weight, TTree* data,
TEventList* training,
TEventList* test,
TNeuron::ENeuronType type = TNeuron::kSigmoid,
const char* extF = "", const char* extD = "");
virtual ~TMultiLayerPerceptron();
void SetData(TTree*);
void SetTrainingDataSet(TEventList* train);
void SetTestDataSet(TEventList* test);
void SetTrainingDataSet(const char* train);
void SetTestDataSet(const char* test);
void SetLearningMethod(TMultiLayerPerceptron::ELearningMethod method);
void SetEventWeight(const char*);
void Train(Int_t nEpoch, Option_t* option = "text");
Double_t Result(Int_t event, Int_t index = 0) const;
Double_t GetError(Int_t event) const;
Double_t GetError(TMultiLayerPerceptron::EDataSet set) const;
void ComputeDEDw() const;
void Randomize() const;
void SetEta(Double_t eta);
void SetEpsilon(Double_t eps);
void SetDelta(Double_t delta);
void SetEtaDecay(Double_t ed);
void SetTau(Double_t tau);
void SetReset(Int_t reset);
inline Double_t GetEta() const { return fEta; }
inline Double_t GetEpsilon() const { return fEpsilon; }
inline Double_t GetDelta() const { return fDelta; }
inline Double_t GetEtaDecay() const { return fEtaDecay; }
inline Double_t GetTau() const { return fTau; }
inline Int_t GetReset() const { return fReset; }
inline TString GetStructure() const { return fStructure; }
inline TNeuron::ENeuronType GetType() const { return fType; }
void DrawResult(Int_t index = 0, Option_t* option = "test") const;
void DumpWeights(Option_t* filename = "-") const;
void LoadWeights(Option_t* filename = "");
Double_t Evaluate(Int_t index, Double_t* params) const;
void Export(Option_t* filename = "NNfunction", Option_t* language = "C++") const;
virtual void Draw(Option_t *option="");
protected:
void AttachData();
void BuildNetwork();
void GetEntry(Int_t) const;
void MLP_Stochastic(Double_t*);
void MLP_Batch(Double_t*);
Bool_t LineSearch(Double_t*, Double_t*);
void SteepestDir(Double_t*);
void ConjugateGradientsDir(Double_t*, Double_t);
void SetGammaDelta(TMatrixD&, TMatrixD&, Double_t*);
bool GetBFGSH(TMatrixD&, TMatrixD &, TMatrixD&);
void BFGSDir(TMatrixD&, Double_t*);
Double_t DerivDir(Double_t*);
Double_t GetCrossEntropyBinary() const;
Double_t GetCrossEntropy() const;
Double_t GetSumSquareError() const;
private:
TMultiLayerPerceptron(const TMultiLayerPerceptron&);
TMultiLayerPerceptron& operator=(const TMultiLayerPerceptron&);
void ExpandStructure();
void BuildFirstLayer(TString&);
void BuildHiddenLayers(TString&);
void BuildLastLayer(TString&, Int_t);
void Shuffle(Int_t*, Int_t) const;
void MLP_Line(Double_t*, Double_t*, Double_t);
TTree* fData;
Int_t fCurrentTree;
Double_t fCurrentTreeWeight;
TObjArray fNetwork;
TObjArray fFirstLayer;
TObjArray fLastLayer;
TObjArray fSynapses;
TString fStructure;
TString fWeight;
TNeuron::ENeuronType fType;
TNeuron::ENeuronType fOutType;
TString fextF;
TString fextD;
TEventList *fTraining;
TEventList *fTest;
ELearningMethod fLearningMethod;
TTreeFormula* fEventWeight;
TTreeFormulaManager* fManager;
Double_t fEta;
Double_t fEpsilon;
Double_t fDelta;
Double_t fEtaDecay;
Double_t fTau;
Double_t fLastAlpha;
Int_t fReset;
Bool_t fTrainingOwner;
Bool_t fTestOwner;
ClassDef(TMultiLayerPerceptron, 4)
};
#endif
Last change: Wed Jun 25 08:49:41 2008
Last generated: 2008-06-25 08:49
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.