{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "89ebe308",
   "metadata": {},
   "source": [
    "# rf901_numintconfig\n",
    "Numeric algorithm tuning: configuration and customization of how numeric (partial) integrals are executed\n",
    "\n",
    "\n",
    "\n",
    "\n",
    "**Author:**  Clemens Lange, Wouter Verkerke (C++ version)  \n",
    "<i><small>This notebook tutorial was automatically generated with <a href= \"https://github.com/root-project/root/blob/master/documentation/doxygen/converttonotebook.py\">ROOTBOOK-izer</a> from the macro found in the ROOT repository  on Tuesday, May 19, 2026 at 08:35 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "bf13b815",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:30.222974Z",
     "iopub.status.busy": "2026-05-19T20:35:30.222833Z",
     "iopub.status.idle": "2026-05-19T20:35:31.192859Z",
     "shell.execute_reply": "2026-05-19T20:35:31.192277Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "65d6c9b3",
   "metadata": {},
   "source": [
    "Adjust global 1D integration precision\n",
    "----------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "99d76224",
   "metadata": {},
   "source": [
    "Print current global default configuration for numeric integration\n",
    "strategies"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "83daa6d3",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:31.196707Z",
     "iopub.status.busy": "2026-05-19T20:35:31.196552Z",
     "iopub.status.idle": "2026-05-19T20:35:31.364167Z",
     "shell.execute_reply": "2026-05-19T20:35:31.363441Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "Requested precision: 1e-07 absolute, 1e-07 relative\n",
      "\n",
      "1-D integration method: RooIntegrator1D (RooImproperIntegrator1D if open-ended)\n",
      "2-D integration method: RooAdaptiveIntegratorND (N/A if open-ended)\n",
      "N-D integration method: RooAdaptiveIntegratorND (N/A if open-ended)\n",
      "\n",
      "Available integration methods:\n",
      "\n",
      "*** RooBinIntegrator ***\n",
      "Capabilities: [1-D] [2-D] [N-D] \n",
      "Configuration: \n",
      "  1)  numBins = 100\n",
      "\n",
      "*** RooIntegrator1D ***\n",
      "Capabilities: [1-D] \n",
      "Configuration: \n",
      "  1)        sumRule = Trapezoid(idx = 0)\n",
      "\n",
      "  2)  extrapolation = Wynn-Epsilon(idx = 1)\n",
      "\n",
      "  3)       maxSteps = 20\n",
      "  4)       minSteps = 999\n",
      "  5)       fixSteps = 0\n",
      "\n",
      "*** RooIntegrator2D ***\n",
      "Capabilities: [2-D] \n",
      "Configuration: \n",
      "(Depends on 'RooIntegrator1D')\n",
      "\n",
      "*** RooSegmentedIntegrator1D ***\n",
      "Capabilities: [1-D] \n",
      "Configuration: \n",
      "  1)  numSeg = 3\n",
      "(Depends on 'RooIntegrator1D')\n",
      "\n",
      "*** RooSegmentedIntegrator2D ***\n",
      "Capabilities: [2-D] \n",
      "Configuration: \n",
      "(Depends on 'RooSegmentedIntegrator1D')\n",
      "\n",
      "*** RooImproperIntegrator1D ***\n",
      "Capabilities: [1-D] [OpenEnded] \n",
      "Configuration: \n",
      "(Depends on 'RooIntegrator1D')\n",
      "\n",
      "*** RooMCIntegrator ***\n",
      "Capabilities: [1-D] [2-D] [N-D] \n",
      "Configuration: \n",
      "  1)   samplingMode = Importance(idx = 0)\n",
      "\n",
      "  2)        genType = QuasiRandom(idx = 0)\n",
      "\n",
      "  3)        verbose = false(idx = 0)\n",
      "\n",
      "  4)          alpha = 1.5\n",
      "  5)    nRefineIter = 5\n",
      "  6)  nRefinePerDim = 1000\n",
      "  7)     nIntPerDim = 5000\n",
      "\n",
      "*** RooAdaptiveIntegratorND ***\n",
      "Capabilities: [2-D] [N-D] \n",
      "Configuration: \n",
      "  1)  maxEval2D = 100000\n",
      "  2)  maxEval3D = 1e+06\n",
      "  3)  maxEvalND = 1e+07\n",
      "  4)    maxWarn = 5\n",
      "\n",
      "*** RooAdaptiveGaussKronrodIntegrator1D ***\n",
      "Capabilities: [1-D] [OpenEnded] \n",
      "Configuration: \n",
      "  1)  maxSeg = 100\n",
      "  2)  method = 21Points(idx = 2)\n",
      "\n",
      "\n",
      "*** RooGaussKronrodIntegrator1D ***\n",
      "Capabilities: [1-D] [OpenEnded] \n",
      "Configuration: \n",
      "\n"
     ]
    }
   ],
   "source": [
    "ROOT.RooAbsReal.defaultIntegratorConfig().Print(\"v\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "5d5d0943",
   "metadata": {},
   "source": [
    "Example: Change global precision for 1D integrals from 1e-7 to 1e-6\n",
    "\n",
    "The relative epsilon (change as fraction of current best integral estimate) and\n",
    "absolute epsilon (absolute change w.r.t last best integral estimate) can be specified\n",
    "separately. For most pdf integrals the relative change criterium is the most important,\n",
    "however for certain non-pdf functions that integrate out to zero a separate absolute\n",
    "change criterium is necessary to declare convergence of the integral\n",
    "\n",
    "NB: ROOT.This change is for illustration only. In general the precision should be at least 1e-7\n",
    "for normalization integrals for MINUIT to succeed."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "f7e979a8",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:31.365939Z",
     "iopub.status.busy": "2026-05-19T20:35:31.365802Z",
     "iopub.status.idle": "2026-05-19T20:35:31.474781Z",
     "shell.execute_reply": "2026-05-19T20:35:31.474039Z"
    }
   },
   "outputs": [],
   "source": [
    "ROOT.RooAbsReal.defaultIntegratorConfig().setEpsAbs(1e-6)\n",
    "ROOT.RooAbsReal.defaultIntegratorConfig().setEpsRel(1e-6)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "813e2b16",
   "metadata": {},
   "source": [
    "Numeric integration of landau pdf\n",
    "------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "a87f89fe",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:31.476614Z",
     "iopub.status.busy": "2026-05-19T20:35:31.476467Z",
     "iopub.status.idle": "2026-05-19T20:35:31.609140Z",
     "shell.execute_reply": "2026-05-19T20:35:31.608329Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", -10, 10)\n",
    "landau = ROOT.RooLandau(\"landau\", \"landau\", x, 0.0, 0.1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8cd6e48a",
   "metadata": {},
   "source": [
    "Disable analytic integration from demonstration purposes"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "da017532",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:31.610958Z",
     "iopub.status.busy": "2026-05-19T20:35:31.610819Z",
     "iopub.status.idle": "2026-05-19T20:35:31.718835Z",
     "shell.execute_reply": "2026-05-19T20:35:31.718127Z"
    }
   },
   "outputs": [],
   "source": [
    "landau.forceNumInt(True)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "db54a9eb",
   "metadata": {},
   "source": [
    "Activate debug-level messages for topic integration to be able to follow\n",
    "actions below"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "890c6a6b",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:31.720876Z",
     "iopub.status.busy": "2026-05-19T20:35:31.720742Z",
     "iopub.status.idle": "2026-05-19T20:35:31.867171Z",
     "shell.execute_reply": "2026-05-19T20:35:31.866447Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/plain": [
       "3"
      ]
     },
     "execution_count": 6,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "ROOT.RooMsgService.instance().addStream(ROOT.RooFit.DEBUG, Topic=ROOT.RooFit.Integration)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "e5901fa0",
   "metadata": {},
   "source": [
    "Calculate integral over landau with default choice of numeric integrator"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "5206e886",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:31.869021Z",
     "iopub.status.busy": "2026-05-19T20:35:31.868885Z",
     "iopub.status.idle": "2026-05-19T20:35:32.002144Z",
     "shell.execute_reply": "2026-05-19T20:35:32.001729Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#3] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>\n",
      "[#3] DEBUG:Integration -- landau: Adding observable x as shape dependent\n",
      "[#3] DEBUG:Integration -- landau: Adding parameter 0 as value dependent\n",
      "[#3] DEBUG:Integration -- landau: Adding parameter 0.1 as value dependent\n",
      "[#3] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)\n",
      "[#3] INFO:Integration -- landau: Observables (x) are numerically integrated\n"
     ]
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      " [1] int_dx landau(x) =  0.09896533620544186\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n"
     ]
    }
   ],
   "source": [
    "intLandau = landau.createIntegral({x})\n",
    "val = intLandau.getVal()\n",
    "print(\" [1] int_dx landau(x) = \", val)  # setprecision(15)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "be6af7e6",
   "metadata": {},
   "source": [
    "Same with custom configuration\n",
    "-----------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a1017417",
   "metadata": {},
   "source": [
    "Construct a custom configuration which uses the adaptive Gauss-Kronrod technique\n",
    "for closed 1D integrals"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "c76c49df",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:32.010243Z",
     "iopub.status.busy": "2026-05-19T20:35:32.010108Z",
     "iopub.status.idle": "2026-05-19T20:35:32.130955Z",
     "shell.execute_reply": "2026-05-19T20:35:32.130490Z"
    }
   },
   "outputs": [],
   "source": [
    "customConfig = ROOT.RooNumIntConfig(ROOT.RooAbsReal.defaultIntegratorConfig())\n",
    "integratorGKNotExisting = customConfig.method1D().setLabel(\"RooAdaptiveGaussKronrodIntegrator1D\")\n",
    "if integratorGKNotExisting:\n",
    "    print(\"WARNING: RooAdaptiveGaussKronrodIntegrator is not existing because ROOT is built without Mathmore support\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "134d9b45",
   "metadata": {},
   "source": [
    "Calculate integral over landau with custom integral specification"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "063a935b",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:32.142265Z",
     "iopub.status.busy": "2026-05-19T20:35:32.142135Z",
     "iopub.status.idle": "2026-05-19T20:35:32.272305Z",
     "shell.execute_reply": "2026-05-19T20:35:32.271571Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#3] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>\n",
      "[#3] DEBUG:Integration -- landau: Adding observable x as shape dependent\n",
      "[#3] DEBUG:Integration -- landau: Adding parameter 0 as value dependent\n",
      "[#3] DEBUG:Integration -- landau: Adding parameter 0.1 as value dependent\n",
      "[#3] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)\n",
      "[#3] INFO:Integration -- landau: Observables (x) are numerically integrated\n"
     ]
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      " [2] int_dx landau(x) =  0.09895710292189495\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooAdaptiveGaussKronrodIntegrator1D to calculate Int(x)\n"
     ]
    }
   ],
   "source": [
    "intLandau2 = landau.createIntegral({x}, NumIntConfig=customConfig)\n",
    "val2 = intLandau2.getVal()\n",
    "print(\" [2] int_dx landau(x) = \", val2)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "ee4204cc",
   "metadata": {},
   "source": [
    "Adjusting default config for a specific pdf\n",
    "-------------------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "063f5602",
   "metadata": {},
   "source": [
    "Another possibility: associate custom numeric integration configuration\n",
    "as default for object 'landau'"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "0dab2dd2",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:32.273679Z",
     "iopub.status.busy": "2026-05-19T20:35:32.273541Z",
     "iopub.status.idle": "2026-05-19T20:35:32.381506Z",
     "shell.execute_reply": "2026-05-19T20:35:32.380612Z"
    }
   },
   "outputs": [],
   "source": [
    "landau.setIntegratorConfig(customConfig)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "581a1edc",
   "metadata": {},
   "source": [
    "Calculate integral over landau custom numeric integrator specified as\n",
    "object default"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "99516d3a",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:32.382974Z",
     "iopub.status.busy": "2026-05-19T20:35:32.382836Z",
     "iopub.status.idle": "2026-05-19T20:35:32.488192Z",
     "shell.execute_reply": "2026-05-19T20:35:32.487458Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      " [3] int_dx landau(x) =  0.09895710292189495\n",
      "[#3] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>\n",
      "[#3] DEBUG:Integration -- landau: Adding observable x as shape dependent\n",
      "[#3] DEBUG:Integration -- landau: Adding parameter 0 as value dependent\n",
      "[#3] DEBUG:Integration -- landau: Adding parameter 0.1 as value dependent\n",
      "[#3] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)\n",
      "[#3] INFO:Integration -- landau: Observables (x) are numerically integrated\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooAdaptiveGaussKronrodIntegrator1D to calculate Int(x)\n"
     ]
    }
   ],
   "source": [
    "intLandau3 = landau.createIntegral({x})\n",
    "val3 = intLandau3.getVal()\n",
    "print(\" [3] int_dx landau(x) = \", val3)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "c1925e36",
   "metadata": {},
   "source": [
    "Another possibility: Change global default for 1D numeric integration\n",
    "strategy on finite domains"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "3898612e",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:35:32.489536Z",
     "iopub.status.busy": "2026-05-19T20:35:32.489415Z",
     "iopub.status.idle": "2026-05-19T20:35:32.609727Z",
     "shell.execute_reply": "2026-05-19T20:35:32.608980Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "--- RooAbsArg ---\n",
      "  Value State: clean\n",
      "  Shape State: clean\n",
      "  Attributes:  [SnapShot_ExtRefClone] \n",
      "  Address: 0x563bb45ef740\n",
      "  Clients: \n",
      "  Servers: \n",
      "  Proxies: \n",
      "--- RooAbsCategory ---\n",
      "  Value = 1 \"15Points)\n",
      "  Possible states:\n",
      "    15Points\t1\n",
      "    21Points\t2\n",
      "    31Points\t3\n",
      "    41Points\t4\n",
      "    51Points\t5\n",
      "    61Points\t6\n",
      "    WynnEpsilon\t0\n"
     ]
    }
   ],
   "source": [
    "if not integratorGKNotExisting:\n",
    "    ROOT.RooAbsReal.defaultIntegratorConfig().method1D().setLabel(\"RooAdaptiveGaussKronrodIntegrator1D\")\n",
    "\n",
    "    # Adjusting parameters of a specific technique\n",
    "    # ---------------------------------------------------------------------------------------\n",
    "\n",
    "    # Adjust maximum number of steps of ROOT.RooIntegrator1D in the global\n",
    "    # default configuration\n",
    "    ROOT.RooAbsReal.defaultIntegratorConfig().getConfigSection(\"RooIntegrator1D\").setRealValue(\"maxSteps\", 30)\n",
    "\n",
    "    # Example of how to change the parameters of a numeric integrator\n",
    "    # (Each config section is a ROOT.RooArgSet with ROOT.RooRealVars holding real-valued parameters\n",
    "    #  and ROOT.RooCategories holding parameters with a finite set of options)\n",
    "    customConfig.getConfigSection(\"RooAdaptiveGaussKronrodIntegrator1D\").setRealValue(\"maxSeg\", 50)\n",
    "    customConfig.getConfigSection(\"RooAdaptiveGaussKronrodIntegrator1D\").setCatLabel(\"method\", \"15Points\")\n",
    "\n",
    "    # Example of how to print set of possible values for \"method\" category\n",
    "    customConfig.getConfigSection(\"RooAdaptiveGaussKronrodIntegrator1D\").find(\"method\").Print(\"v\")"
   ]
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3",
   "language": "python",
   "name": "python3"
  },
  "language_info": {
   "codemirror_mode": {
    "name": "ipython",
    "version": 3
   },
   "file_extension": ".py",
   "mimetype": "text/x-python",
   "name": "python",
   "nbconvert_exporter": "python",
   "pygments_lexer": "ipython3",
   "version": "3.12.12"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 5
}
