Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
tmva101_Training.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_tmva
3## \notebook -nodraw
4## This tutorial show how you can train a machine learning model with any package
5## reading the training data directly from ROOT files. Using XGBoost, we illustrate
6## how you can convert an externally trained model in a format serializable and readable
7## with the fast tree inference engine offered by TMVA.
8##
9## \macro_code
10## \macro_output
11##
12## \date August 2019
13## \author Stefan Wunsch
14
15# XGBoost has to be imported before ROOT to avoid crashes because of clashing
16# std::regexp symbols that are exported by cppyy.
17# See also: https://github.com/wlav/cppyy/issues/227
18from xgboost import XGBClassifier
19
20import ROOT
21import numpy as np
22
23from tmva100_DataPreparation import variables
24
25
26def load_data(signal_filename, background_filename):
27 # Read data from ROOT files
28 data_sig = ROOT.RDataFrame("Events", signal_filename).AsNumpy()
29 data_bkg = ROOT.RDataFrame("Events", background_filename).AsNumpy()
30
31 # Convert inputs to format readable by machine learning tools
32 x_sig = np.vstack([data_sig[var] for var in variables]).T
33 x_bkg = np.vstack([data_bkg[var] for var in variables]).T
34 x = np.vstack([x_sig, x_bkg])
35
36 # Create labels
37 num_sig = x_sig.shape[0]
38 num_bkg = x_bkg.shape[0]
39 y = np.hstack([np.ones(num_sig), np.zeros(num_bkg)])
40
41 # Compute weights balancing both classes
42 num_all = num_sig + num_bkg
43 w = np.hstack([np.ones(num_sig) * num_all / num_sig, np.ones(num_bkg) * num_all / num_bkg])
44
45 return x, y, w
46
47if __name__ == "__main__":
48 # Load data
49 x, y, w = load_data("train_signal.root", "train_background.root")
50
51 # Fit xgboost model
52 bdt = XGBClassifier(max_depth=3, n_estimators=500)
53 bdt.fit(x, y, sample_weight=w)
54
55 # Save model in TMVA format
56 print("Training done on ",x.shape[0],"events. Saving model in tmva101.root")
57 ROOT.TMVA.Experimental.SaveXGBoost(bdt, "myBDT", "tmva101.root", num_inputs=x.shape[1])
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...