#ifndef ROOT_TMVA_MethodBDT
#define ROOT_TMVA_MethodBDT
#include <vector>
#ifndef ROOT_TH2
#include "TH2.h"
#endif
#ifndef ROOT_TTree
#include "TTree.h"
#endif
#ifndef ROOT_TMVA_MethodBase
#include "TMVA/MethodBase.h"
#endif
#ifndef ROOT_TMVA_DecisionTree
#include "TMVA/DecisionTree.h"
#endif
#ifndef ROOT_TMVA_Event
#include "TMVA/Event.h"
#endif
namespace TMVA {
   class SeparationBase;
   class MethodBDT : public MethodBase {
   public:
      
      MethodBDT( const TString& jobName,
                 const TString& methodTitle,
                 DataSetInfo& theData,
                 const TString& theOption = "",
                 TDirectory* theTargetDir = 0 );
      
      MethodBDT( DataSetInfo& theData,
                 const TString& theWeightFile,
                 TDirectory* theTargetDir = NULL );
      virtual ~MethodBDT( void );
      virtual Bool_t HasAnalysisType( Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets );
      
      
      void InitEventSample();
      
      virtual std::map<TString,Double_t> OptimizeTuningParameters(TString fomType="ROCIntegral", TString fitType="FitGA");
      virtual void SetTuneParameters(std::map<TString,Double_t> tuneParameters);
      
      void Train( void );
      
      void Reset( void );
      using MethodBase::ReadWeightsFromStream;
      
      void AddWeightsXMLTo( void* parent ) const;
      
      void ReadWeightsFromStream( std::istream& istr );
      void ReadWeightsFromXML(void* parent);
      
      void WriteMonitoringHistosToFile( void ) const;
      
      Double_t GetMvaValue( Double_t* err = 0, Double_t* errUpper = 0);
      
      UInt_t   GetNTrees() const {return fForest.size();}
   private:
      Double_t GetMvaValue( Double_t* err, Double_t* errUpper, UInt_t useNTrees );
      Double_t PrivateGetMvaValue( const TMVA::Event *ev, Double_t* err=0, Double_t* errUpper=0, UInt_t useNTrees=0 );
      void     BoostMonitor(Int_t iTree);
   public:
      const std::vector<Float_t>& GetMulticlassValues();
      
      const std::vector<Float_t>& GetRegressionValues();
      
      Double_t Boost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
      
      const Ranking* CreateRanking();
      
      void DeclareOptions();
      void ProcessOptions();
      void SetMaxDepth(Int_t d){fMaxDepth = d;}
      void SetMinNodeSize(Double_t sizeInPercent);
      void SetMinNodeSize(TString sizeInPercent);
      void SetNTrees(Int_t d){fNTrees = d;}
      void SetAdaBoostBeta(Double_t b){fAdaBoostBeta = b;}
      void SetNodePurityLimit(Double_t l){fNodePurityLimit = l;} 
      void SetShrinkage(Double_t s){fShrinkage = s;}
      void SetUseNvars(Int_t n){fUseNvars = n;}
      void SetBaggedSampleFraction(Double_t f){fBaggedSampleFraction = f;}
      
      inline const std::vector<TMVA::DecisionTree*> & GetForest() const;
      
      inline const std::vector<const TMVA::Event*> & GetTrainingEvents() const;
      inline const std::vector<double> & GetBoostWeights() const;
      
      std::vector<Double_t> GetVariableImportance();
      Double_t GetVariableImportance(UInt_t ivar);
      Double_t TestTreeQuality( DecisionTree *dt );
      
      void MakeClassSpecific( std::ostream&, const TString& ) const;
      
      void MakeClassSpecificHeader( std::ostream&, const TString& ) const;
      void MakeClassInstantiateNode( DecisionTreeNode *n, std::ostream& fout,
                                     const TString& className ) const;
      void GetHelpMessage() const;
   protected:
      void DeclareCompatibilityOptions();
   private:
      
      void Init( void );
      void PreProcessNegativeEventWeights();
      
      Double_t AdaBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
      
      Double_t AdaCost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
      
      Double_t Bagging( );
      
      Double_t RegBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt );
      
      Double_t AdaBoostR2( std::vector<const TMVA::Event*>&, DecisionTree *dt );
      
      
      
      Double_t GradBoost( std::vector<const TMVA::Event*>&, DecisionTree *dt, UInt_t cls = 0);
      Double_t GradBoostRegression(std::vector<const TMVA::Event*>&, DecisionTree *dt );
      void InitGradBoost( std::vector<const TMVA::Event*>&);
      void UpdateTargets( std::vector<const TMVA::Event*>&, UInt_t cls = 0);
      void UpdateTargetsRegression( std::vector<const TMVA::Event*>&,Bool_t first=kFALSE);
      Double_t GetGradBoostMVA(const TMVA::Event *e, UInt_t nTrees);
      void     GetBaggedSubSample(std::vector<const TMVA::Event*>&);
      Double_t GetWeightedQuantile(std::vector<std::pair<Double_t, Double_t> > vec, const Double_t quantile, const Double_t SumOfWeights = 0.0);
      std::vector<const TMVA::Event*>       fEventSample;     
      std::vector<const TMVA::Event*>       fValidationSample;
      std::vector<const TMVA::Event*>       fSubSample;       
      std::vector<const TMVA::Event*>      *fTrainSample;     
      Int_t                           fNTrees;          
      std::vector<DecisionTree*>      fForest;          
      std::vector<double>             fBoostWeights;    
      Double_t                        fSigToBkgFraction;
      TString                         fBoostType;       
      Double_t                        fAdaBoostBeta;    
      TString                         fAdaBoostR2Loss;  
      Double_t                        fTransitionPoint; 
      Double_t                        fShrinkage;       
      Bool_t                          fBaggedBoost;     
      Bool_t                          fBaggedGradBoost; 
      Double_t                        fSumOfWeights;    
      std::map< const TMVA::Event*, std::pair<Double_t, Double_t> >       fWeightedResiduals;  
      std::map< const TMVA::Event*,std::vector<double> > fResiduals; 
      
      SeparationBase                 *fSepType;         
      TString                         fSepTypeS;        
      Int_t                           fMinNodeEvents;   
      Float_t                         fMinNodeSize;     
      TString                         fMinNodeSizeS;    
      Int_t                           fNCuts;           
      Bool_t                          fUseFisherCuts;   
      Double_t                        fMinLinCorrForFisher; 
      Bool_t                          fUseExclusiveVars; 
      Bool_t                          fUseYesNoLeaf;    
      Double_t                        fNodePurityLimit; 
      UInt_t                          fNNodesMax;       
      UInt_t                          fMaxDepth;        
      DecisionTree::EPruneMethod       fPruneMethod;     
      TString                          fPruneMethodS;    
      Double_t                         fPruneStrength;   
      Double_t                         fFValidationEvents;    
      Bool_t                           fAutomatic;       
      Bool_t                           fRandomisedTrees; 
      UInt_t                           fUseNvars;        
      Bool_t                           fUsePoissonNvars; 
      UInt_t                           fUseNTrainEvents; 
      Double_t                         fBaggedSampleFraction;     
      TString                          fNegWeightTreatment;     
      Bool_t                           fNoNegWeightsInTraining; 
      Bool_t                           fInverseBoostNegWeights; 
      Bool_t                           fPairNegWeightsGlobal;   
      Bool_t                           fTrainWithNegWeights; 
      Bool_t                           fDoBoostMonitor; 
      
      TTree*                           fMonitorNtuple;   
      Int_t                            fITree;           
      Double_t                         fBoostWeight;     
      Double_t                         fErrorFraction;   
      Double_t                         fCss;             
      Double_t                         fCts_sb;          
      Double_t                         fCtb_ss;          
      Double_t                         fCbb;             
      
      Bool_t                           fDoPreselection;  
      std::vector<Double_t>            fVariableImportance; 
      void                             DeterminePreselectionCuts(const std::vector<const TMVA::Event*>& eventSample);
      Double_t                         ApplyPreselectionCuts(const Event* ev);
      
      std::vector<Double_t> fLowSigCut;
      std::vector<Double_t> fLowBkgCut;
      std::vector<Double_t> fHighSigCut;
      std::vector<Double_t> fHighBkgCut;
      
      std::vector<Bool_t>  fIsLowSigCut;  
      std::vector<Bool_t>  fIsLowBkgCut;  
      std::vector<Bool_t>  fIsHighSigCut; 
      std::vector<Bool_t>  fIsHighBkgCut; 
      
      Bool_t fHistoricBool; 
      
      static const Int_t               fgDebugLevel;     
      
      ClassDef(MethodBDT,0)  
   };
} 
const std::vector<TMVA::DecisionTree*>& TMVA::MethodBDT::GetForest()         const { return fForest; }
const std::vector<const TMVA::Event*> & TMVA::MethodBDT::GetTrainingEvents() const { return fEventSample; }
const std::vector<double>&              TMVA::MethodBDT::GetBoostWeights()   const { return fBoostWeights; }
#endif