{
"cells": [
{
"cell_type": "markdown",
"id": "0024c455",
"metadata": {},
"source": [
"# tmva004_RStandardScaler\n",
"This tutorial illustrates the usage of the standard scaler as preprocessing\n",
"method.\n",
"\n",
"\n",
"\n",
"\n",
"**Author:** Stefan Wunsch \n",
"This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, March 19, 2024 at 07:20 PM."
]
},
{
"cell_type": "code",
"execution_count": 1,
"id": "d684d746",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2024-03-19T19:20:09.608115Z",
"iopub.status.busy": "2024-03-19T19:20:09.607363Z",
"iopub.status.idle": "2024-03-19T19:20:10.554033Z",
"shell.execute_reply": "2024-03-19T19:20:10.552792Z"
}
},
"outputs": [],
"source": [
"using namespace TMVA::Experimental;"
]
},
{
"cell_type": "markdown",
"id": "3b6fbc8c",
"metadata": {},
"source": [
"Load data used to fit the parameters"
]
},
{
"cell_type": "code",
"execution_count": 2,
"id": "47e5763a",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2024-03-19T19:20:10.562040Z",
"iopub.status.busy": "2024-03-19T19:20:10.561663Z",
"iopub.status.idle": "2024-03-19T19:20:15.810748Z",
"shell.execute_reply": "2024-03-19T19:20:15.809485Z"
}
},
"outputs": [],
"source": [
"ROOT::RDataFrame df(\"TreeS\", \"http://root.cern/files/tmva_class_example.root\");\n",
"auto x = AsTensor(df);"
]
},
{
"cell_type": "markdown",
"id": "1ae4bdd5",
"metadata": {},
"source": [
"Create standard scaler and fit to data"
]
},
{
"cell_type": "code",
"execution_count": 3,
"id": "235bf6d6",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2024-03-19T19:20:15.815860Z",
"iopub.status.busy": "2024-03-19T19:20:15.815552Z",
"iopub.status.idle": "2024-03-19T19:20:16.238722Z",
"shell.execute_reply": "2024-03-19T19:20:16.237273Z"
}
},
"outputs": [],
"source": [
"RStandardScaler scaler;\n",
"scaler.Fit(x);"
]
},
{
"cell_type": "markdown",
"id": "15a6b1e4",
"metadata": {},
"source": [
"Compute transformation"
]
},
{
"cell_type": "code",
"execution_count": 4,
"id": "25cfc1ff",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2024-03-19T19:20:16.254945Z",
"iopub.status.busy": "2024-03-19T19:20:16.254502Z",
"iopub.status.idle": "2024-03-19T19:20:16.714088Z",
"shell.execute_reply": "2024-03-19T19:20:16.698726Z"
}
},
"outputs": [],
"source": [
"auto y = scaler.Compute(x);"
]
},
{
"cell_type": "markdown",
"id": "a6fd535d",
"metadata": {},
"source": [
"Plot first variable scaled and unscaled"
]
},
{
"cell_type": "code",
"execution_count": 5,
"id": "153f1b45",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2024-03-19T19:20:16.722118Z",
"iopub.status.busy": "2024-03-19T19:20:16.721721Z",
"iopub.status.idle": "2024-03-19T19:20:17.926644Z",
"shell.execute_reply": "2024-03-19T19:20:17.925150Z"
}
},
"outputs": [
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAAxgAAALoCAIAAAB0+sSeAAAABmJLR0QAAAAAAAD5Q7t/AAAgAElEQVR4nO3dwZKzOIIuUHyjHqs6wOuZB+pegVczD9S9BsfMe3EX6tRQgElbNsayzomKivydMimM0/pSEtJpHMcKAIDH/b+jKwAAkCtBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQKI/nj/EMAxN0zx/nDc4nU5HVwEAvt84jkdX4U1SglTXdZfLJbxGMZr0fZ9LnNpQzoUHgJ0U1W3x8NDeMAyXy6Wu66qquq6rqqrv+7quz+fzyyu3h3HT0bUDAHKSEqSm/6/rumma8E8AgKKkBKnQHVVV1fV6nQ7niVMAQFEeDlJN01yv12EYwrje9P9fMEcKAOB+p4SJQXESWV3XIVGFWVOf3yN1OqWcLwBwv6Ja28RTDZkpdEFNv/5wRV1aADhEUa3tw6d6KzZ1XRcG+D5ZUZcWAA5RVGubMtl8GZjCmgivqREAQCYeWJAzTDMPXxe11hYAwKrHOt+GH8tOKXOkAICqsNY25VQz2lxvpqhLCwCHKKq1fequvZnPT1e/jkiWc+EBYCdFBamUTYtvxZEsXrUsKgkAZOHhu/biRsV2/AUACvdwkKqqKmxU/OqaAABkJmWvvR2qAQCQn8QgJU4BADw82XwYhrAs53LKuWlSAEBRHg5STdP0fb9HVQAA8lLQSg9VYStbAMAhimptU+7aq6qq67qmacLoXtg05pWVAgDIQUqQaprmcrmEr0OKOp/Py933AAC+28NBKkw27/s+9kJ1Xde2bYxWAACFeHgUMyx8EFLU6XTq+z48Mv36Y9lrDwD2Zo7Ulo2o9OEpKljubGOXGwAgTUqQul6vsxlRWUQoAIDXSul867puOSPq88f1qsI6G+E7/DYg/wI+FeC1impt00+167owU6ppmlxu2Svq0sIXeEOKqgQpeLWiWtuCTrUq7NLCFwhBar/f2r2PD2UqqrV9eI5U13Wn06lpGotwAgCFSwlSbdtWVXU+nyUqAKBkKSubh9lR4zhOE1Uu06QAAF7lBaOY8Sa+zx8QLWrUFr6AOVKQo6Ja2z+SnzkMQ9d11+u1qqq6rj9/7QMAgNdK3LT4dDqdz+fr9dq27TiOIVS9um4A8ErDMGxM7e267nPasucrE072NbXhtpTJ5jE/jeP4Oe+5O502HV07APZ1vV5vBanL5bJcbvool8vlyXu5hmEIo0bs6uGhvY8K7AnKGbUFAPZ2b4/UMAzTaLz8p+4cAKA09/ZIhV6oGJ5m/wSAL9N1XdM0YRu0W1uiTacIr26Ydudz44/bqM/wI+FnsaPxPnVd13V96599399/qHiEWIe6rvu+Xy0WlqqKbhXr+356wDiFa+bRSgK/qqrd/9u78pQjtFYbbcS0mQglYzsVv4gFpu3O8rvhZ91qc1e/O63Y7Giz1nDZgq9W5qhWr6jW9oAgtfpuWA1Jy7fF6i/A6ttxtT5FXVp4gzekKEGqIPu/dR4NUstwEwuEQ02bwtWnx6YttH2x/LKdWj49Fo7NXPzu7Gizf46TBnT1TPdWVGubsvzBM4ZhCPdETN+a4Xqfz+dpybhI1exduLyRITxx+YZz2ye8xzju+x8caDpGFpqh0AbFQbT43b7vY09BeFbbtrFAeCTeRlfX9axbYbXvYPrccfLL0HVdXdfxaKFhnc1dvuv0eN6dgetVPVLhjTJ97vTxZfCf/d2wfHp8I06LLcN7POY9lQTulHuPTu7151GP9kjNWqvQ4oQ/72NDs3q0R4fVls1Wtdl9NavMRsN6fx1eqKjW9t09UiE+r06Rq34y9fLx2T+nC2PE/q1psdlfAACQYGNko2maOE4SViJ8tMXpui6uYjgbk1k1W/gwtpjLvrFfK88LPbCO1PV6na1xkLzkwa/X+9bbMRabLdi6LB/6PPVtAjDzqqYhLKwYN0wL63n2fX9PggkNaNhgLdpeP3N1QvC0WZx9Swv4HvcGqXCNX/VTl+vWr17vjQHj7WX+q6r69R0JQGl+bcg22p2l0BLFxqjrusvlEtcgWArfiqsetG37UCfWsvLTH7Rs8jSC7/HYOlLPCx1Fy/fZcnWN7SPc+eO8jQCYCeMVs1wS/nl/l0H4c33a/xSCVPz6fD7P/uCfTV+Z/ayNBqtt23Cj1fQpoU9r/JkgNTsp3VFv88c//vGPR5/zX//1X8k/L7y3QpaaLmIW30DT98EzfWBN06xumZQ8HDm6dwggf33fn8/n8/kch9Xi7eTVI70GoTk7n8/h1rx4kHCEpmlCuIlrY4ZZUGF4LrRQsR2cVWB1GvHlcln+rOlNgquV4UM9Ob99dR2p2YoXt27uG/96f9+tu/Nufev5ygNTud/1lnv9SfPrYpjB8sHpjXLjWnO2eqf56ndnzw3fCl/fuhFv42jLkworXR/V6hXV2p7+/PPP1bQ087//+7/x6/Hpvpk4NS9c7LAu/rSXMqTvuq6XnZOhWBhaHoYhZPxllWIX6/Rbp9Pp+coDUejhzfe3Kvf686TYxDwzAPLrQW7dVVetzRi+88dtPCvhmC9XVGv7x//8z//8Wugf//hHDFL//Oc/n/+py3nij47m3jmY/dC0QQCK8pLA8etBttdQeOePYw+/rCP1r3/963Q6/fd//3dVVX/++ec4jv/xH/+xRz1CkIq5J7wPNmbe/XrT36+ZHQDgSVtB6m9/+9t//ud/hq//+c9/3tN39auw/tgy30zn6FW3F8ZYXS9q9fEQwgQpAGA/60EqdESF4by///3vL+yIil1N04S0eptemIg3W+x1uY55XOv81wMCALzWynSwv/3tb9MZUS8fy4uLHYTOpDh+t1wNNi5VEJfhD/+c1Xl6wDD7Ku52vDxgOdPf4A1yn6yde/3hMxXV2v7lVP/xj3+E6VBVVf39739/Zr2obctlx2+tqT8ruXofX/XXZdA2DljUpYU3yD2I5F5/+ExFtbb/PtV//etfcTrUn3/++ZLpUL+6/77TO2eO/3rAoi4tvEHuQST3+sNnKqq1PY3jOO2I2mMs73MUdWnhDXIPIrnXHz5TUa3tH3Ee0ts6ogAAvkPKxnP5xsyiMjK8Qe49OrnXHz5TUa3tH0dX4N22Ny0u58IDAM/74yVbvmREVAIAXqWgzreqsM5GeIPch8Zyrz98pqJa21/22gMA4BZBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJ77f1FOSuxAu+w+YHzGj614FDFBSlRCXiTN6Qo4GjFBSmAt9rvjzdBDT6AOVIAAIkEKQCARIIUAEAiQQoAIJEgBQCQSJACAEgkSAEAJBKkAAASCVIAAImKW9ncXnsAwKsUF6REJQDgVYoLUlAge7IB7MQcKfhyUhTAfvRIQRGMaQPsQZACSrdTp53sCiUwtAcAkEiPFFCufUc8zU6DAuiRAgBIpEcKKJvbGoEn6JECCiZFAc/RIwUUb5+pUiGkuXcPvltxQcpeewDAqxQXpEQlAOBVzJECAEgkSAEAJDosSDVNc/rRNM0wDKvFuq47TdwqNgzD9IBd1+1WcQCAfzu9f87QMAzn83n5eNu2swDUNM31ev212K0DLk/tdDrgfOFY/753LN83/htWKNj1rr39XvnsLy1fq6jW9oAeqRB66roef7RtW1XV5XKZdjh1XRdSVN/3G8WWB+z7PjzeNM1bTgjYjXWegM/27swYe49mPzd0Pk17m8I6BbP+p1CsruuYpbquu1wuswPe+ilFZWQI8u62yLn2eqQoVlGt7bt7pEIAqut69njoPVpOgZqN4oV/Tsf7QooKnVWzoy2fDgDwQsdMNl/OfJpFqFsBKCakX8uHrHZrcjoAJbvzhqdnjv/kMbuu0xeQhXcHqfi2mE5gitOhZm+aZcdVFN6gG29TE6QAWAo3g4dGp67ruq6v1+v5fH5hqzEMw/V6fTJIDcMQhlz4cAesbD6OY3gTz3Zr6fv+VofTVHjT3/mz7i8JQAlCOpm2ONXPBFydQCQ4YGjv/pD+zN8Ht557SpVcEwA+RMhJbdvO2ohfRznglgMmm99a/uB8Ps/exHu8p8dUL68JAIdY/Ut7NV11Xdc0TdM0y56qsBB09GuDNS2/WnjjZ/HRkoNFmjDtaZqignjb3fSfy2LjT6Bp23acLBm1LDY7YHzuC84BslJVY8Zv/Jxrv3vdc35xjrLaNGyXjJZ//8/EVQ9D2xTaqWA55XfWwH1IG/0q+dY8wbt7pFYnlU8fubMXKvzd8OvY38Z0dQA+yun04v+WpksVbvT9xInesbEM03NjCzX7bkhOt44WbqiaJqdwtNndV9PglfT6cZB3pLWJjR8avhUS/UZX0+zx6bOmQoSavi/HwjIyBHl3W+Rcez1Sjwon9ML/btnubRp/WpBZyxKL9X1f1/Xyu7GJmfVIrTZn0weXBTYawSzkW/MEx6wj9euyBY+uF7V8PHR9WQQBIBfj+OL/bum6bhzHvu/btg2NSLiRPLY4qy3IOI6hQJjkNFvEZ/vUVof2whfhmLeWlSYDbw5ut8anl3Onbk1yqv7azxRj+/SPg/gGnf2U958vHC7vbouca69HKiOzFqe6MUl3anXqSPjWtEcqNlK3muDwo5fjKoe00a+Sb80THLAbTlxKoK7rkOvjak+zysSScbvi1WJh/Y/lAWfLhFSF7f4DQd4bsuVce3vtfaBb25FVP01JaDhC63OrvZjefh5vxJs+JRQIe8WGr+u6Xu21Cm1WLDz91nYdPlxZre0h8e3X8eloFvlvFVsecJnux8IyMgR5d1vkXHs9Uh9oo+Gb9gytFmvbNnQyhYZp48jLOVKrjdd0TvCywIFt9PPyrXmCI0+1/3FnyecPWNSlhSDv1jbn2gtSH+jWEjzjX4NLCFXLaSThicsgNZsbPgtSt4JXrMYyM92/TMNnyrfmCQo61bGwSwtB3q1tzrUXpD5THLho2zb87R1Ty3IBglAmFgh/qMeVDmdPr/5643k8WoxZq0cbJ7EpFJgOxbzzlXmhfGueoKBTHQu7tBDk3drmXHtB6mOtzhOfLZezLLZ6S1P81rRTarkg53LK+WzwZDZBJQwj5ttm5VvzBCVNB5vMXr+lqFeDQuQ9Iznn2pts/uHilPPttQZmKx3c+a2Eo8UqfcHaB0VNNi/oVKvCLi0Eebe2OddekKJYRbW2xyzICQDwBQQpAIBEghQAQKI/jq4AAE/47R6aZxUz0wXS6JECAEikRwogT3v3Fe3d1wVfQY8UAEAiQQoAIJEgBQCQSJACAEhU3GTz7e32ylnSHgB4XnFBSlQCAF7F0B4AQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgUXFbxNhrDwB4leKClKgEALyKoT0AgESCFABAIkEKACCRIAUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACBRcXvt2bQYAHiV4oKUqAQAvIqhPQCARIIUAECidw/tNU2zXaDrummZrusul0v8Z9/3q0cYhqHruuv1Gv7Ztm3Xdc/WFQBg0+nNc4a253pXf81ATdPEbLRaIBiG4Xw+Lw+1PLXT6d3nC4cLv3O5vvFzrn3Oda+q6gtOgMMU1dq++1Q3OopCz1Psc4p9URuPBCGc1XU9DEM1yVXxkWnJci4tBHm3hjnXPue6V1X1BSfAYYpqbT/lVENImvY2hXg0638KfVTThBTT1fREYpaanV1RlxaCvFvDnGufc92rqvqCE+AwRbW2H3GqIffMOpBCkJpVb5mQVvPWrceLurQQ5N0a5lz7nOteVdUXnACHKaq1/Yi79kI2mqaoWyOAcURvNma3LF/X9bIYAMALHR+kQgZq23b5rRCGVsXpULcK/Hp7IADAk44PUmGG0/JGvFvlN9LV0vKmPwCAVzk4SIV+o77vN777zJGXTqmSawIAfKuD99oLPUa3Qs8eM5zKmf4GAOztyB6pjdlRG31Rd47WmWYOAOztyCC1OjvqTiFp/Tr299CEKgCAhxwWpEKP0a2gExLSRufTLEIt+5/CI+7dAwD2c1iQCh1Rt4LOo+tFLR/fnn0FAPC8w9YeXV24fGp175fleuVxrfPpBnxxt2NbxEDey1PnXPuc615V1RecAIcpqrU95lRv7YU3ExcdCBPSQ65aPivGprqum6YZhiH8c7a3cVXYpYUg79Yw59rnXPeqqr7gBDhMUa3tMacaeptmm+utiiEpuPWU2H0VLVNUVdilhSDv1jDn2udc96qqvuAEOExRrW0ep3rnzPGYsTZW48zifOGF8m4Nc659znWvquoLToDDFNXaFnSqVWGXFoK8W8Oca59z3auq+oIT4DBFtbbH77UHAJCpg7eIeb/tXfPKSdAAwPOKC1KiEgDwKob2AAASCVIAAIkEKQCARIIUAEAiQQoAIJEgBQCQSJACAEgkSAEAJBKkAAASCVIAAImK2yLGXnsAwKsUF6REJQDgVQztAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIVt0WMvfYAgFcpLkiJSgDAqxjaAwBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkKi4BTmBF9vcLQDgu+mRAp4gRQFlK65Hyl578Hp+cYBSFRekRCUA4FUM7QEAJBKkAAASHRakhmE4TQzDsFqs67p7ig3D0DRNLNZ13W4VBwD4t9Mhc4aaprler7MH67qe5aTVYm3bznLSMAzn83n5U5andjodc75woHB/xV5v/H2PnrfsX5vsT4DDFNXaHtAj1XVdiEd934/jOI5j3/dVVV2v12mQWhZr27aqqsvlMstbIUXVdT09WlVVTdO85YQAgEIdkBnDAgR930+DTtd1l8tl2ikVis36n0If1bRYeGL11/6n2Ec1O7uiMjIEeqSOkv1rk/0JcJiiWtt3n+qtiFNVVQhMMTaFIDUrtnz6at669XhRlxYCQeoo2b822Z8AhymqtX33OlIh1oRButVvLb+eip1YYXb5Rvm6rmdjhQAAr/XuIBWmPYUMNAxDCDpN06zOZ6rr+tZxQpDayEmrE9UBAF7osJXNp1u1hElO01lTGwkpdDXd+VNkKQBgP8esIxXmObVt2/d927ah5+l8Pi+XP0j+Ebeee0qVXBMA4Fsd1iMV+5/i/6/X6/l8nt189/KfW870NwBgb8f0SLVtO+sxur8v6s7ROtPMAYC9vTtIbcwfD+4MQNPerGd+HABAsg/dtDgkpI3Op+0OrfiIxc0BgP28O0iFBZ9u5Z5q0dU0K7m6XtTq49N1FgAA9vDuIBW7mmbRJ+6XFx8Ji3bOdiMOCyVM1/MMx5mtvRnzkyAFAOzngEXc4+541U9yikN4y63xwhdxu+LVYnHtzbquwyqdcbfjWZAqatF6CGwRc5TsX5vsT4DDFNXaHnOq0ywVTPchnpotUH6r2PKAyxRVFXZpIRCkjpL9a5P9CXCYolrbI091Ni/q15J3FtsoWdSlhUCQOkr2r032J8BhimptCzrVqrBLC4EgdZTsX5vsT4DDFNXafujyBwAAn0+QAgBIdNhee0fZ3n64nK5IPopNsQEyVVyQEpX4NFIUQL6KC1LwmSR8gByZIwUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkKi4lc3ttQcAvEpxQUpUAgBexdAeAEAiQQoAIJEgBQCQSJACAEgkSAEAJBKkAAASCVIAAIkEKQCARIIUAECi4lY2B3inzV2pXsBmDXCs4oKUvfYAgFcpLkiJSsB77P1hs3dfF3APc6QAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiYrbIsZeewAP2HUnGh+55K+4ICUqAQCvUlyQAuAuu/7ZactlvoU5UgAAiQ4IUsMwNDcMwzAr3HXdaWJZYHrMWKzrup1PAgCgOr1/zlDXdZfLZfVbbdtOM1DTNNfrdbtMVVXDMJzP5+XRlqd2Oh1wvrAtDHHk+sbMu/Z5y/u1z7v2/KKo1vaAOVKhV6mu66ZpZt+aPtJ1XUhRfd+Hx0MCu1wuofsqlgwpqq7rcOSYq1a7uAAAXuWAzBgWIIjxaLvYah9VzEzVpH9reiIxS83OrqiMTC7y/ss879rnLe/XPu/a84uiWtvDgtSvP3e12DIhreatW48XdWnJRd4NSt61z1ver33etecXRbW2H3rX3q3Z4rETazZmtyxf1/WyGADAC707SMVkM73J7tZ9diEMbRxnIydtjxsCADzvsCAVJpLHqHS5XKabt2wkpI10tbS86Q8A4FUOG9rr+34cx2EYxnHs+z48OOuXeqZX6dZzT6mSawIAfKt3B6mu60KEmgadpmnatq2qara+1B4znMZUL68JAJC7A3qkVvuK7u+LunO0zjRzAGBvH3fX3p0BKCStX8f+HppQBQDwkAOG9rY3wpsmpI3Op1mEWsav8Ih79wCA/bw7SIU9Xm7lnujR9aKWj4cQJkgBAPt5d5AKk8qXewyHR8J3N0qG2ejTYiFCXa/XaeSK+UmQAgD2c9gWMdXPvsXDMMQ1pWb9T7Hk7J6+WZ3DBnzLAy638ytq0XpykfdWGXnXPm95v/Z5155fFNXaHnOqMfpEy83yVksuw1YQty6OVjdFLurSkou8G5S8a5+3vF/7vGvPL4pqbY881el2MfeUvLPYRsmiLi25yLtBybv2ecv7tc+79vyiqNa2oFOtCru05CLvBiXv2uct79c+79rzi6Ja249bRwoAIBd/HF2Bd9veNa+cBA0APK+4ICUqAQCvYmgPACCRIAUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkKi4LWLstQcAvEpxQUpUAgBexdAeAEAiQQoAIJEgBQCQSJACAEgkSAEAJBKkAAASCVIAAIkEKQCARIIUAEAiQQoAIJEgBQCQqLi99mxaDAC8SnFBSlQCAF7F0B4AQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgUXFbxNhrDwB4leKClKgEALyKoT0AgESCFABAoo8IUqfT6dbUpa7rThPDMKwWG4ahaZpYrOu6/WoLABCcDp8z1DTN9Xqt1mYvxW9NtW07y0nDMJzP5+WRlwc8nY4/X5gJf0Tk+sbMu/Z5y/u1z7v2/KKo1vbgHqmu65ZRafatvu/HcRzHsW3bqqoul8usXyqkqLquQ7G+78PjTdPsWXcAoHRHZsZZT9KsJmGwb9b/FPqo6rqOWarrusvlMnt6PPLymOVkZHKR91/medc+b3m/9nnXnl8U1doeeaohKvV9fyv0LB9cJqTVvHXr8aIuLbnIu0HJu/Z5y/u1z7v2/KKo1vawob0w7ta27eoA3K3Z4rHwbHRvWb6u62UxAIAXOiZIDcMQRui2b68LYejWEarNnGSCFACwt2OCVBie24hBG9/aSFdLt2ayAwA874AgFfqK4r11v5Z85qcsnVIl1wQA+Fbv3msvLGpQ1/U9IWmPGU7lTH8DAPb21h6pYRjCUgW/JqSNmHXnaJ1p5gDA3o6ZI7U6avbQ7i4haf3arfXQhCoAgId8xF57SyEhbXQ+zSLUsv8pPOLePQBgP28NUk3TjGvCd8PXoUfq0fWilo+HECZIAQD7+dAeqaqqws56s92IwxSr8K0gRKjr9TqNXDE/CVIAwH4+YhH31d1g4uPVT3IKKWpZMmzAV1VVuBkwrPZZVVXf97MgVdSi9eQi760y8q593vJ+7fOuPb8oqrX9iFO9FaSqSUgKptsVT8Wti6NliqoKu7TkIu8GJe/a5y3v1z7v2vOLolrbPE71zpnjMWNtrMaZxflSlLwblLxrn7e8X/u8a88vimptCzrVqrBLSy7yblDyrn3e8n7t8649vyiqtf3cyeYAAB/u3VvEHG5717xyEjQA8LzigpSoBPAp9t4P3gc++zO0BwCQqLgeKQCOt3df0d59XfBDjxQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBExS3Iaa89AOBVigtSohLwTexWB8cytAcAkKi4HikokX3HvpHd6uAT6JGCb6c9BNiNHikog6kuADvQIwUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgETFLX9grz0SWIkJgFXFBSlRiUdJUQDcUlyQgjQSOABL5kgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkEiQAgBIJEgBACQqbmVze+0BAK9SXJASlQCAVzG0BwCQSJACAEgkSAEAJDomSA3D0DTN6UfTNMMwrJbsuu40cavY7IBd1+1WdwCAfzu9f/J113WXy2X5eNu2swDUNM31ev212DAM5/N5ecDlqZ1OB5wvuQs3emb8xsn+BDhG3m+cvGufvaJa2wN6pEKKatt2/NG2bXw86roupKi+72fFZv1SIUXVdR2K9X0fHm+a5g2nAwCUa3yvEIZi6IlCZabpavnIOI51XS+5gG0AABEoSURBVM+eHg44O5GYpZY/5flToDRVNeb9xsn+BDhG3m+cvGufvaJa23f3SIXOpOUcppCQllOgZiXDP6fjfbF/a1os9kWZLAUA7OfdQappmrZtbw26/RqAYoFZ5Lo/mQEAvMq7VzZfTUjDMIROplnACmFoVbhNbyMnrU5UBwB4oSPXkYpLG4TZ4n3f3+pwmtpIV0uyFACwnyOD1K/jbs/cdnfruadUyTUBAL7VwUFqHMe+78NU8fP5PItWe8xwSp6W//KaAAC5O36LmKZpuq6LWSo+eKv8naN1ppkDAHs7PkgFj65TEJLWr2N/D02oAgB4yLuD1J3zjUJC2uh8mkWoZf9TeMTi5gDAfg7btHj7kUfXi1o+vrqeAgDAC707SIXcs9xjODwyXaB8NmsqWK5jHtc6n0aumJ8EKQBgPwfszxyH9uq6DotqxiG8WWViydmuxrNice3N2QGnC1PFA7r/jkdlv4t89ifAMfJ+4+Rd++wV1doec6rLZcfbtl2dbz4rWdf16u14XdfFmBUsU1RV2KXlVbL/QM7+BDhG3m+cvGufvaJa2yNPNUaiXwfg7pw5/usBi7q0vEr2H8jZnwDHyPuNk3fts1dUa1vQqVaFXVpeJfsP5OxPgGPk/cbJu/bZK6q1/ZR1pAAAsvPH0RV4t+1VrMpJ0ADA84oLUqISAPAqhvYAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiYrbIsZee3yizbclAB+ruCAlKvFxpCiAbBUXpOBDifgAGTJHCgAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgUXErm9trDwB4leKClKgEALyKoT0AgESCFABAIkEKACCRIAUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAouK2iLHXHgDwKsUFKVEJAHgVQ3sAAIkEKQCARIcFqaZpTj+aphmGYbVY13WniVvFhmGYHrDrut0qDgDwb6f3zxkahuF8Pi8fr+t6lpOaprler7NibdvOctKtAy5P7XQ64HzJXbg/Ycc3zu4/AFLk/cbMu/bZK6q1PaBHKoSeuq7HH33fV1V1vV6nCanrupCi+r4Pxdq2rarqcrnM8tbsgOFoVVU1TfOWEwLgI51O+/4H7++R6rrucrlUi+6i2KsUHw/rFMz6n0If1bTvavWAy6PFY5aTkXkVPVKUKe835ntSTq6vzu6Kam3f3SMVAlDoW5qKvUez3qbZKF7453S8L6So2QHj0UyWAijROO77H/w4ZrL5r4NutwLQnXmrqqq6rpfFAABe6IAeqXEcl0EqJp7pt0IYunWcajMnmSAFAOztU9aRClOa4gjdRkLaSFdLy5v++FamkwLwfscHqbBSVFVVdV3PRuie6VW69dxTquSa8AauDwCHOHKvven6T8vVoap9ZjiVcx9BgVxbAN7syJXN4/pPfd/f3xd152idaeYAz7MME2w7pkcqjpT1fZ82fhee9etzH5pQBQDwkAN6pOKMqNXb94Lw+Ebn0+yJy/6n8Ih79wDSWIYJ7vHuIBWG8Jbb6s08ul7U8vEQwgQpAGA/7w5SYSHye/JNWAphthvxch3zuNb5NHLF4wtSAMB+3r0bzq/rCEzrEwvH7YqXZaqfDfiqqqrrummaYRjibsezIFXU7j9F+YZNwXKtPaTI/l2f/Qnsq6jW9qODVDUJScGtMcG4dXG0Oo29qEtblLw/0/KuPaTI/l2f/Qnsq6jWNo9TvXPm+Oo+M1NFXdqi5P2ZlnftIUX27/rsT2BfRbW2BZ1qVdilLUren2l51x5SZP+uz/4E9lVUa3v8FjEAAJk6couYQ2xP0ionQQMAzysuSIlKAMCrGNoDAEgkSAEAJBKkAAASCVIAAIkEKQCARIIUAEAiQQoAIJEgBQCQSJACAEgkSAEAJCpuixh77QEAr1JckBKVAIBXMbQHAJBIkAIASCRIAQAkEqQAABIJUgAAiQQpAIBEghQAQCJBCgAgkSAFAJBIkAIASCRIAQAkKm6vPZsWAwCvUlyQEpUAgFcxtAcAkEiQAgBIJEgBACQSpAAAEglSAACJirtrDwBeY3M9nWe5xzwTeqQAABLpkQKAB+3aXbRrRxevJkjBHXyuAbDG0B78RooC4IbieqTstUci7w0AFooLUqISAPAqhvYAABIJUgAAiQ4OUsMwnE6nYRhuFei67jRxq+QwDE3TxGJd1+1TXwCA/3M6ds5Q0zTX67Xv+6Zpbn139mDbtrOcNAzD+XxePn15aqfTwefLTsItBHtd232PDiXyW7Ul/1enqNb2yB6pruuWOWn53b7vx3Ecx7Ft26qqLpfLrF8qpKi6rkOxvu/D46vhDADgVQ4IUnG07nK5bBQL323bNuahruvqug5fTI8WvojpqmmakKU2UhoAwPM+fbL5bBQv/HOakGLemhabZq9dqwcAlOyYHqnxx0aZ1cdjQpqN7i3Lh76rjWnsAABP+ugeqRCGVoWEtJGTTJACAPb2oUFqIyFtpKsl06QAgP18aJAKnulVuvXcU6rkmgAA3+qj99rbY4ZTOStbAAB7+9AeqY2+qDtH60wzBwD29qFB6lchaf069vfQhCoAgId8aJAKCWmj82kWoZb9T+ER9+4BAPv56CBV3b1e1PLxEMIEKQBgPx8apKqfxcpnuxEv1zGPa51PI1fMT4LU5ziddvwPyNGuHws+GXiPg/dnDssK9H2/mnjiogNxu+Lwz1mdm6YJ/U91XTdNMwxD3O14dtii9qP+KO/5RNvr2ua/Ezt8mrw/E/aW/2dOUa3tRwepahKSgrquV2/H67putgXy6jGLurQfJe+PhbxrDyXK+7c279pXVWGtbR6neufM8ZixNlbjzOJ8v0/eHwt51x5KlPdvbd61r6rCWtuCTrUq7NJ+lLw/FvKuPZQo79/avGtfVYW1tp872RwA4MN99BYxe9jeNa+cBA0APK+4ICUqAQCvYmgPACCRIAUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkKi4LWLstQcAvEpxQUpUAgBexdAeAEAiQQoAIJEgBQCQqLg5UnytzdsIAGAPeqT4ClIUAEfQI8UXcUsmAO+lRwoAIJEeKQD4PHvPWNCF/yJ6pAAAEumRAoBPsndfkbtzXqq4IGWvPQDgVYoLUqISQCHMMuINigtS3KKvFwAeJUhRVVIU8F3MMuJtBCn+j25qAHiI5Q8AABIJUgAAiQSpr7W90ANZK+riOtlvVdTJVuWdbzkEKQCARIIUAEAiQQoAIJHlDwCgPHvO2SpqLZ3igpS99gCAVykuSIlKABTNuu8vVVyQ4jCF/WoBX+/RT7WHyvurPxcmm/MWUhQA30iP1JbT6fTQUOCnlX/Q6aEJgimV2fvF2fXgn1T+UQ8dv6iT3bt8USebUP4hH1X5cUw4/r0fgadTtfcH8kdd2dx9SY9U0zSnH03THF2dHZ1O9/73aGEA4FHf0CM1uxHver2eTqe+778sUe0dd8YHf8aj5QHg+2TfIxXTUt/34ziO41jXdVVV5/P5yGrtZhzv/e+hwmMlEgHAw/IOUsMwXK/XqqrGcYyJahiG8EXXdcdUK18PJK/xlFAeAL5L9kGqqqrQBTXVtm01SVQAAHv4hiC1nAsVHgmdVW92//zuhPng4yNPeKjw+18oADbs15r4yH+tvINUiEq3gtSb7T8Z3HsfAD7LN9y1d79db0z7WWHj0bjzYPldVwN6rCoAvN44VvcvOhU8VD72YD3k/vKlNSVlBakvsL3p8jOFSyv/UZX5tPIfVZncy39UZT6t/EdV5tPK712Z3f/sL0ZhQSr/hVmzPwEA+CJ5z5ECADhQ3kEqLHywXObAwgcAwBvkHaQCQQoAOETeQSqsXb5cLyoEqbAsJwDATh67Q/4DhfsU2raNG8IMwxA22sv91ACAD5d3j1RVVX3fV1V1uVyapum6rmmakKIe6o46nU6P3ziajWEYmqY5/Wia5uuHPks732gYhtPp9B3nO7uIR1fnfb7pIq4q6tezwI/f6Lsb1r8Y8xey1FTbtvc/PW7Vt1sFj3QrUD70EmVk+WYI6ro+umrvEN7Mfd8fXZFnrV7ELzive3zNRVwq7deztI/fqe9uWGe+5yT7Hw89a/pG36deB1v+3sZTPq5SO1p+LsfP7q//8IpXNvc2OH4ExxMp50P5ay7iqtJ+PUv7+I2+vmGdKeIkb5n9eXR0dV4vvJuXf+196yfXrQ+peKEPqdXeln/1Zt0G37pY3/qmDb7sIq4q7deztI/f6Osb1qXs50g9I8ymutXb/AXCYHychh/dWn8rd7fu1owzbL7vlL9PuEaxCyoKl9UVzFdpv56lffxGX9+wLpUbpMJvb9u2XzyPtWmajRP81hP/1vO6peu6+IfR0XV5gdDALC9ieGS51sl3+LKLuKGcX8+SP36/u2FdcUQ32PFCWI6drkW9FPEPhe8bO7jlW8cOlr7gym6cgov4lcr59Ry//eO32Ia10B6p0Pf4xZ2rS13XhZtRY79rOX8xJKyIAbxHCb+ehXz8FtiwBiUGqfAOLmoEtyryzV39fH5VVVXX9XKyAnCgcn49S/j4LbNhDYoLUl3XXa/Xuq6/8g+CDcMwjOPY9334y+98Pn/373ZY1fByuVRV1bbtd58s5KW0X8+v//gttmH9t6PHFl9geYPPrTPdvrP6iLo/7P6T3ZDRWiZp5xufVdd1dnMRnrnE4VvZnfLUxink8qZ90hdcxG1Z/3o+L6OP3zt9QcP6pD/u+cj+cAl9wqvr1i+37ftAL6lb13Xhb8HP98zFzXQiwie//d6gruvr9Ro21pg+/mV/wRcr91/P52X08fuofBvWJ31DkCrqt/Ghkw3v4DHnG6ofvbhxykW+7W5R7+dblpcv3wtK9AW/nvf7go9f7vQNQep+TdOsvq2/+x1fzh/34S+eQj6mv1XXdefzeble1K3lHMlFmb+eJXz8ltmwThU32bwoYS5CuCV16lvvNw4d5np0shYv33QgYBiGEK2+eHTg65X261nax2/JTiWkxV99cXCOg9bhforYIFXfeL6rI/RT33fKM+EVyH32yTAMobGZvWm/e5pF9B0XcanAX8+iPn6XvrhhnRGkqurbr3fTNLOBkm9tkAr8pJ75mjY4ZqnoW9+0S19zEWfK/PUs5+N36bsb1ilBqhRxYP7LPp35Yt60fAfv5O8mSAEAJDLZHAAgkSAFAJBIkAIASCRIAQAkEqQAABIJUgDAXNM0q+s1NE1TyFJYdxKkAIAV1+t1lqWWS4xiHSkAYMVsnf2u68KeiZLDlCAFAKyImzWFqFDOri8P+ePoCgAAn6hpmrqupwN8bdseWqNPpEcKALgpbjhd13XcN5DIZHMA4Ka+78MXbtZbpUcKALhJj9Q2PVIAwLrQC1XXdVVV1+tVkFrSIwUArJjetTe7g2+1cFVVq2t4fjdBCgBYMVtHKqzGeWuA73Q6lTn2Z2gPAJiLg3qxkymEpNUBvgI7oiJBCgD4i2EYwiLms8wU7uALY3xR13Ul7xtjaA8ASBTmTvV9fz6fDe0BADzgfD63bWtoDwDgMWEPmcIX6rTXHgDwsGEYrterCUJ6pACAh4XpUKcfVVVdr9fT6VTaNCmTzQGAFNPMFCabd11X2nwpQ3sAQIplZiotRVWG9gAAkhnaAwBIpEcKACCRIAUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkEiQAgBIJEgBACQSpAAAEglSAACJBCkAgESCFABAIkEKACCRIAUAkEiQAgBIJEgBACQSpAAAEv1/mSwkNSb2/aMAAAAASUVORK5CYII=\n",
"text/plain": [
""
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"TH1F h1(\"h1\", \";x_{4};N_{Events}\", 20, -4, 4);\n",
"TH1F h2(\"h2\", \";x_{4};N_{Events}\", 20, -4, 4);\n",
"for (std::size_t i = 0; i < x.GetShape()[0]; i++) {\n",
" h1.Fill(x(i, 3));\n",
" h2.Fill(y(i, 3));\n",
"}\n",
"h1.SetLineWidth(2);\n",
"h1.SetLineColor(kRed);\n",
"h2.SetLineWidth(2);\n",
"h2.SetLineColor(kBlue);\n",
"\n",
"gStyle->SetOptStat(0);\n",
"auto c = new TCanvas(\"\", \"\", 800, 800);\n",
"h2.Draw(\"HIST\");\n",
"h1.Draw(\"HIST SAME\");\n",
"\n",
"TLegend legend(0.7, 0.7, 0.89, 0.89);\n",
"legend.SetBorderSize(0);\n",
"legend.AddEntry(\"h1\", \"Unscaled\", \"l\");\n",
"legend.AddEntry(\"h2\", \"Scaled\", \"l\");\n",
"legend.Draw();\n",
"\n",
"c->DrawClone();"
]
}
],
"metadata": {
"kernelspec": {
"display_name": "ROOT C++",
"language": "c++",
"name": "root"
},
"language_info": {
"codemirror_mode": "text/x-c++src",
"file_extension": ".C",
"mimetype": " text/x-c++src",
"name": "c++"
}
},
"nbformat": 4,
"nbformat_minor": 5
}