Sophie

Sophie

distrib > Fedora > 13 > i386 > media > updates > by-pkgid > f3eb4c16ba6256fe5a10e54bf649f01f > files > 1236

python-mpmath-doc-0.17-1.fc13.noarch.rpm

<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN"
  "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">

<html xmlns="http://www.w3.org/1999/xhtml">
  <head>
    <meta http-equiv="Content-Type" content="text/html; charset=utf-8" />
    
    <title>Number-theoretical, combinatorial and integer functions &mdash; mpmath v0.17 documentation</title>
    <link rel="stylesheet" href="../_static/default.css" type="text/css" />
    <link rel="stylesheet" href="../_static/pygments.css" type="text/css" />
    <script type="text/javascript">
      var DOCUMENTATION_OPTIONS = {
        URL_ROOT:    '../',
        VERSION:     '0.17',
        COLLAPSE_MODINDEX: false,
        FILE_SUFFIX: '.html',
        HAS_SOURCE:  true
      };
    </script>
    <script type="text/javascript" src="../_static/jquery.js"></script>
    <script type="text/javascript" src="../_static/doctools.js"></script>
    <link rel="top" title="mpmath v0.17 documentation" href="../index.html" />
    <link rel="up" title="Mathematical functions" href="index.html" />
    <link rel="next" title="q-functions" href="qfunctions.html" />
    <link rel="prev" title="Zeta functions, L-series and polylogarithms" href="zeta.html" /> 
  </head>
  <body>
    <div class="related">
      <h3>Navigation</h3>
      <ul>
        <li class="right" style="margin-right: 10px">
          <a href="../genindex.html" title="General Index"
             accesskey="I">index</a></li>
        <li class="right" >
          <a href="../modindex.html" title="Global Module Index"
             accesskey="M">modules</a> |</li>
        <li class="right" >
          <a href="qfunctions.html" title="q-functions"
             accesskey="N">next</a> |</li>
        <li class="right" >
          <a href="zeta.html" title="Zeta functions, L-series and polylogarithms"
             accesskey="P">previous</a> |</li>
        <li><a href="../index.html">mpmath v0.17 documentation</a> &raquo;</li>
          <li><a href="index.html" accesskey="U">Mathematical functions</a> &raquo;</li> 
      </ul>
    </div>  

    <div class="document">
      <div class="documentwrapper">
        <div class="bodywrapper">
          <div class="body">
            
  <div class="section" id="number-theoretical-combinatorial-and-integer-functions">
<h1>Number-theoretical, combinatorial and integer functions<a class="headerlink" href="#number-theoretical-combinatorial-and-integer-functions" title="Permalink to this headline">¶</a></h1>
<p>For factorial-type functions, including binomial coefficients,
double factorials, etc., see the separate
section <a class="reference external" href="gamma.html"><em>Factorials and gamma functions</em></a>.</p>
<div class="section" id="fibonacci-numbers">
<h2>Fibonacci numbers<a class="headerlink" href="#fibonacci-numbers" title="Permalink to this headline">¶</a></h2>
<div class="section" id="fibonacci-fib">
<h3><tt class="xref docutils literal"><span class="pre">fibonacci()</span></tt>/<tt class="xref docutils literal"><span class="pre">fib()</span></tt><a class="headerlink" href="#fibonacci-fib" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.fibonacci">
<tt class="descclassname">mpmath.</tt><tt class="descname">fibonacci</tt><big>(</big><em>n</em>, <em>**kwargs</em><big>)</big><a class="headerlink" href="#mpmath.fibonacci" title="Permalink to this definition">¶</a></dt>
<dd><p><tt class="docutils literal"><span class="pre">fibonacci(n)</span></tt> computes the <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>-th Fibonacci number, <img class="math" src="../_images/math/2663041fd2b9fbd85c85ad832f67c7a743a09c52.png" alt="F(n)"/>. The
Fibonacci numbers are defined by the recurrence <img class="math" src="../_images/math/3c80df4dced24d2e083ecaf5eba04730e2abd010.png" alt="F(n) = F(n-1) + F(n-2)"/>
with the initial values <img class="math" src="../_images/math/e1676ed7f6f03f442bbb3620cc13f1d86a018db0.png" alt="F(0) = 0"/>, <img class="math" src="../_images/math/0faeea7ef7f0f0e89604b913937e23bd424c7035.png" alt="F(1) = 1"/>. <a title="mpmath.fibonacci" class="reference internal" href="#mpmath.fibonacci"><tt class="xref docutils literal"><span class="pre">fibonacci()</span></tt></a>
extends this definition to arbitrary real and complex arguments
using the formula</p>
<div class="math">
<p><img src="../_images/math/bf1d4599e53249947baa0f7b49dd66f187d02ac2.png" alt="F(z) = \frac{\phi^z - \cos(\pi z) \phi^{-z}}{\sqrt 5}" /></p>
</div><p>where <img class="math" src="../_images/math/2c175f60eecef1de7560c3bdea495d69f26f719d.png" alt="\phi"/> is the golden ratio. <a title="mpmath.fibonacci" class="reference internal" href="#mpmath.fibonacci"><tt class="xref docutils literal"><span class="pre">fibonacci()</span></tt></a> also uses this
continuous formula to compute <img class="math" src="../_images/math/2663041fd2b9fbd85c85ad832f67c7a743a09c52.png" alt="F(n)"/> for extremely large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>, where
calculating the exact integer would be wasteful.</p>
<p>For convenience, <tt class="xref docutils literal"><span class="pre">fib()</span></tt> is available as an alias for
<a title="mpmath.fibonacci" class="reference internal" href="#mpmath.fibonacci"><tt class="xref docutils literal"><span class="pre">fibonacci()</span></tt></a>.</p>
<p><strong>Basic examples</strong></p>
<p>Some small Fibonacci numbers are:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">10</span><span class="p">):</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="n">fibonacci</span><span class="p">(</span><span class="n">i</span><span class="p">))</span>
<span class="gp">...</span>
<span class="go">0.0</span>
<span class="go">1.0</span>
<span class="go">1.0</span>
<span class="go">2.0</span>
<span class="go">3.0</span>
<span class="go">5.0</span>
<span class="go">8.0</span>
<span class="go">13.0</span>
<span class="go">21.0</span>
<span class="go">34.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fibonacci</span><span class="p">(</span><span class="mi">50</span><span class="p">)</span>
<span class="go">12586269025.0</span>
</pre></div>
</div>
<p>The recurrence for <img class="math" src="../_images/math/2663041fd2b9fbd85c85ad832f67c7a743a09c52.png" alt="F(n)"/> extends backwards to negative <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">10</span><span class="p">):</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="n">fibonacci</span><span class="p">(</span><span class="o">-</span><span class="n">i</span><span class="p">))</span>
<span class="gp">...</span>
<span class="go">0.0</span>
<span class="go">1.0</span>
<span class="go">-1.0</span>
<span class="go">2.0</span>
<span class="go">-3.0</span>
<span class="go">5.0</span>
<span class="go">-8.0</span>
<span class="go">13.0</span>
<span class="go">-21.0</span>
<span class="go">34.0</span>
</pre></div>
</div>
<p>Large Fibonacci numbers will be computed approximately unless
the precision is set high enough:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">fib</span><span class="p">(</span><span class="mi">200</span><span class="p">)</span>
<span class="go">2.8057117299251e+41</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">45</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fib</span><span class="p">(</span><span class="mi">200</span><span class="p">)</span>
<span class="go">280571172992510140037611932413038677189525.0</span>
</pre></div>
</div>
<p><a title="mpmath.fibonacci" class="reference internal" href="#mpmath.fibonacci"><tt class="xref docutils literal"><span class="pre">fibonacci()</span></tt></a> can compute approximate Fibonacci numbers
of stupendous size:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fibonacci</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">25</span><span class="p">)</span>
<span class="go">3.49052338550226e+2089876402499787337692720</span>
</pre></div>
</div>
<p><strong>Real and complex arguments</strong></p>
<p>The extended Fibonacci function is an analytic function. The
property <img class="math" src="../_images/math/955ea4f7586d40c1fc76f92e7d448a88fc5b3eff.png" alt="F(z) = F(z-1) + F(z-2)"/> holds for arbitrary <img class="math" src="../_images/math/b13f21416d84e13708696f34dea81026cda583c9.png" alt="z"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fib</span><span class="p">(</span><span class="n">pi</span><span class="p">)</span>
<span class="go">2.1170270579161</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fib</span><span class="p">(</span><span class="n">pi</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span> <span class="o">+</span> <span class="n">fib</span><span class="p">(</span><span class="n">pi</span><span class="o">-</span><span class="mi">2</span><span class="p">)</span>
<span class="go">2.1170270579161</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fib</span><span class="p">(</span><span class="mi">3</span><span class="o">+</span><span class="mi">4</span><span class="n">j</span><span class="p">)</span>
<span class="go">(-5248.51130728372 - 14195.962288353j)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fib</span><span class="p">(</span><span class="mi">2</span><span class="o">+</span><span class="mi">4</span><span class="n">j</span><span class="p">)</span> <span class="o">+</span> <span class="n">fib</span><span class="p">(</span><span class="mi">1</span><span class="o">+</span><span class="mi">4</span><span class="n">j</span><span class="p">)</span>
<span class="go">(-5248.51130728372 - 14195.962288353j)</span>
</pre></div>
</div>
<p>The Fibonacci function has infinitely many roots on the
negative half-real axis. The first root is at 0, the second is
close to -0.18, and then there are infinitely many roots that
asymptotically approach <img class="math" src="../_images/math/7894a934b6faeaabe216e25c212c048c59c74aa7.png" alt="-n+1/2"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">findroot</span><span class="p">(</span><span class="n">fib</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.2</span><span class="p">)</span>
<span class="go">-0.183802359692956</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">findroot</span><span class="p">(</span><span class="n">fib</span><span class="p">,</span> <span class="o">-</span><span class="mi">2</span><span class="p">)</span>
<span class="go">-1.57077646820395</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">findroot</span><span class="p">(</span><span class="n">fib</span><span class="p">,</span> <span class="o">-</span><span class="mi">17</span><span class="p">)</span>
<span class="go">-16.4999999596115</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">findroot</span><span class="p">(</span><span class="n">fib</span><span class="p">,</span> <span class="o">-</span><span class="mi">24</span><span class="p">)</span>
<span class="go">-23.5000000000479</span>
</pre></div>
</div>
<p><strong>Mathematical relationships</strong></p>
<p>For large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>, <img class="math" src="../_images/math/faca141ea81723e491a81891dda2e718140eaa17.png" alt="F(n+1)/F(n)"/> approaches the golden ratio:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">50</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fibonacci</span><span class="p">(</span><span class="mi">101</span><span class="p">)</span><span class="o">/</span><span class="n">fibonacci</span><span class="p">(</span><span class="mi">100</span><span class="p">)</span>
<span class="go">1.6180339887498948482045868343656381177203127439638</span>
<span class="gp">&gt;&gt;&gt; </span><span class="o">+</span><span class="n">phi</span>
<span class="go">1.6180339887498948482045868343656381177203091798058</span>
</pre></div>
</div>
<p>The sum of reciprocal Fibonacci numbers converges to an irrational
number for which no closed form expression is known:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">nsum</span><span class="p">(</span><span class="k">lambda</span> <span class="n">n</span><span class="p">:</span> <span class="mi">1</span><span class="o">/</span><span class="n">fib</span><span class="p">(</span><span class="n">n</span><span class="p">),</span> <span class="p">[</span><span class="mi">1</span><span class="p">,</span> <span class="n">inf</span><span class="p">])</span>
<span class="go">3.35988566624318</span>
</pre></div>
</div>
<p>Amazingly, however, the sum of odd-index reciprocal Fibonacci
numbers can be expressed in terms of a Jacobi theta function:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">nsum</span><span class="p">(</span><span class="k">lambda</span> <span class="n">n</span><span class="p">:</span> <span class="mi">1</span><span class="o">/</span><span class="n">fib</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">n</span><span class="o">+</span><span class="mi">1</span><span class="p">),</span> <span class="p">[</span><span class="mi">0</span><span class="p">,</span> <span class="n">inf</span><span class="p">])</span>
<span class="go">1.82451515740692</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">sqrt</span><span class="p">(</span><span class="mi">5</span><span class="p">)</span><span class="o">*</span><span class="n">jtheta</span><span class="p">(</span><span class="mi">2</span><span class="p">,</span><span class="mi">0</span><span class="p">,(</span><span class="mi">3</span><span class="o">-</span><span class="n">sqrt</span><span class="p">(</span><span class="mi">5</span><span class="p">))</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span><span class="o">**</span><span class="mi">2</span><span class="o">/</span><span class="mi">4</span>
<span class="go">1.82451515740692</span>
</pre></div>
</div>
<p>Some related sums can be done in closed form:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">nsum</span><span class="p">(</span><span class="k">lambda</span> <span class="n">k</span><span class="p">:</span> <span class="mi">1</span><span class="o">/</span><span class="p">(</span><span class="mi">1</span><span class="o">+</span><span class="n">fib</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">k</span><span class="o">+</span><span class="mi">1</span><span class="p">)),</span> <span class="p">[</span><span class="mi">0</span><span class="p">,</span> <span class="n">inf</span><span class="p">])</span>
<span class="go">1.11803398874989</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">phi</span> <span class="o">-</span> <span class="mf">0.5</span>
<span class="go">1.11803398874989</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">f</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">k</span><span class="p">:(</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span><span class="o">**</span><span class="p">(</span><span class="n">k</span><span class="o">+</span><span class="mi">1</span><span class="p">)</span> <span class="o">/</span> <span class="nb">sum</span><span class="p">(</span><span class="n">fib</span><span class="p">(</span><span class="n">n</span><span class="p">)</span><span class="o">**</span><span class="mi">2</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">1</span><span class="p">,</span><span class="nb">int</span><span class="p">(</span><span class="n">k</span><span class="o">+</span><span class="mi">1</span><span class="p">)))</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">nsum</span><span class="p">(</span><span class="n">f</span><span class="p">,</span> <span class="p">[</span><span class="mi">1</span><span class="p">,</span> <span class="n">inf</span><span class="p">])</span>
<span class="go">0.618033988749895</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">phi</span><span class="o">-</span><span class="mi">1</span>
<span class="go">0.618033988749895</span>
</pre></div>
</div>
<p><strong>References</strong></p>
<ol class="arabic simple">
<li><a class="reference external" href="http://mathworld.wolfram.com/FibonacciNumber.html">http://mathworld.wolfram.com/FibonacciNumber.html</a></li>
</ol>
</dd></dl>

</div>
</div>
<div class="section" id="bernoulli-numbers-and-polynomials">
<h2>Bernoulli numbers and polynomials<a class="headerlink" href="#bernoulli-numbers-and-polynomials" title="Permalink to this headline">¶</a></h2>
<div class="section" id="bernoulli">
<h3><tt class="xref docutils literal"><span class="pre">bernoulli()</span></tt><a class="headerlink" href="#bernoulli" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.bernoulli">
<tt class="descclassname">mpmath.</tt><tt class="descname">bernoulli</tt><big>(</big><em>n</em><big>)</big><a class="headerlink" href="#mpmath.bernoulli" title="Permalink to this definition">¶</a></dt>
<dd><p>Computes the nth Bernoulli number, <img class="math" src="../_images/math/9be2b1de79fff66db67de3b18303bdd2912b7ae6.png" alt="B_n"/>, for any integer <img class="math" src="../_images/math/e2c7d62172af0e738dd12f47f3f8ad8ef70cd65f.png" alt="n \ge 0"/>.</p>
<p>The Bernoulli numbers are rational numbers, but this function
returns a floating-point approximation. To obtain an exact
fraction, use <a title="mpmath.bernfrac" class="reference internal" href="#mpmath.bernfrac"><tt class="xref docutils literal"><span class="pre">bernfrac()</span></tt></a> instead.</p>
<p><strong>Examples</strong></p>
<p>Numerical values of the first few Bernoulli numbers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">15</span><span class="p">):</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="s">&quot;</span><span class="si">%s</span><span class="s"> </span><span class="si">%s</span><span class="s">&quot;</span> <span class="o">%</span> <span class="p">(</span><span class="n">n</span><span class="p">,</span> <span class="n">bernoulli</span><span class="p">(</span><span class="n">n</span><span class="p">)))</span>
<span class="gp">...</span>
<span class="go">0 1.0</span>
<span class="go">1 -0.5</span>
<span class="go">2 0.166666666666667</span>
<span class="go">3 0.0</span>
<span class="go">4 -0.0333333333333333</span>
<span class="go">5 0.0</span>
<span class="go">6 0.0238095238095238</span>
<span class="go">7 0.0</span>
<span class="go">8 -0.0333333333333333</span>
<span class="go">9 0.0</span>
<span class="go">10 0.0757575757575758</span>
<span class="go">11 0.0</span>
<span class="go">12 -0.253113553113553</span>
<span class="go">13 0.0</span>
<span class="go">14 1.16666666666667</span>
</pre></div>
</div>
<p>Bernoulli numbers can be approximated with arbitrary precision:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">50</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bernoulli</span><span class="p">(</span><span class="mi">100</span><span class="p">)</span>
<span class="go">-2.8382249570693706959264156336481764738284680928013e+78</span>
</pre></div>
</div>
<p>Arbitrarily large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> are supported:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bernoulli</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">20</span> <span class="o">+</span> <span class="mi">2</span><span class="p">)</span>
<span class="go">3.09136296657021e+1876752564973863312327</span>
</pre></div>
</div>
<p>The Bernoulli numbers are related to the Riemann zeta function
at integer arguments:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="o">-</span><span class="n">bernoulli</span><span class="p">(</span><span class="mi">8</span><span class="p">)</span> <span class="o">*</span> <span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="p">)</span><span class="o">**</span><span class="mi">8</span> <span class="o">/</span> <span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">fac</span><span class="p">(</span><span class="mi">8</span><span class="p">))</span>
<span class="go">1.00407735619794</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">zeta</span><span class="p">(</span><span class="mi">8</span><span class="p">)</span>
<span class="go">1.00407735619794</span>
</pre></div>
</div>
<p><strong>Algorithm</strong></p>
<p>For small <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> (<img class="math" src="../_images/math/67569ba3ff78135cda24dd931ebb5817a602a571.png" alt="n &lt; 3000"/>) <a title="mpmath.bernoulli" class="reference internal" href="#mpmath.bernoulli"><tt class="xref docutils literal"><span class="pre">bernoulli()</span></tt></a> uses a recurrence
formula due to Ramanujan. All results in this range are cached,
so sequential computation of small Bernoulli numbers is
guaranteed to be fast.</p>
<p>For larger <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>, <img class="math" src="../_images/math/9be2b1de79fff66db67de3b18303bdd2912b7ae6.png" alt="B_n"/> is evaluated in terms of the Riemann zeta
function.</p>
</dd></dl>

</div>
<div class="section" id="bernfrac">
<h3><tt class="xref docutils literal"><span class="pre">bernfrac()</span></tt><a class="headerlink" href="#bernfrac" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.bernfrac">
<tt class="descclassname">mpmath.</tt><tt class="descname">bernfrac</tt><big>(</big><em>n</em><big>)</big><a class="headerlink" href="#mpmath.bernfrac" title="Permalink to this definition">¶</a></dt>
<dd><p>Returns a tuple of integers <img class="math" src="../_images/math/e148d95e9a2349eff0d0ccd1ec50ace5fa25269b.png" alt="(p, q)"/> such that <img class="math" src="../_images/math/ed0e40214b4c65391930d730e70ced3aac62cc97.png" alt="p/q = B_n"/> exactly,
where <img class="math" src="../_images/math/9be2b1de79fff66db67de3b18303bdd2912b7ae6.png" alt="B_n"/> denotes the <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>-th Bernoulli number. The fraction is
always reduced to lowest terms. Note that for <img class="math" src="../_images/math/49d98ffe996c0f0d3b82b284cba2113439fefea9.png" alt="n &gt; 1"/> and <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> odd,
<img class="math" src="../_images/math/d72198ae318470c88b2114f13302cfcf4e9e9563.png" alt="B_n = 0"/>, and <img class="math" src="../_images/math/6457ca2dcd40dc21d9aa0db18c66bbac9af4d899.png" alt="(0, 1)"/> is returned.</p>
<p><strong>Examples</strong></p>
<p>The first few Bernoulli numbers are exactly:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">15</span><span class="p">):</span>
<span class="gp">... </span>    <span class="n">p</span><span class="p">,</span> <span class="n">q</span> <span class="o">=</span> <span class="n">bernfrac</span><span class="p">(</span><span class="n">n</span><span class="p">)</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="s">&quot;</span><span class="si">%s</span><span class="s"> </span><span class="si">%s</span><span class="s">/</span><span class="si">%s</span><span class="s">&quot;</span> <span class="o">%</span> <span class="p">(</span><span class="n">n</span><span class="p">,</span> <span class="n">p</span><span class="p">,</span> <span class="n">q</span><span class="p">))</span>
<span class="gp">...</span>
<span class="go">0 1/1</span>
<span class="go">1 -1/2</span>
<span class="go">2 1/6</span>
<span class="go">3 0/1</span>
<span class="go">4 -1/30</span>
<span class="go">5 0/1</span>
<span class="go">6 1/42</span>
<span class="go">7 0/1</span>
<span class="go">8 -1/30</span>
<span class="go">9 0/1</span>
<span class="go">10 5/66</span>
<span class="go">11 0/1</span>
<span class="go">12 -691/2730</span>
<span class="go">13 0/1</span>
<span class="go">14 7/6</span>
</pre></div>
</div>
<p>This function works for arbitrarily large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">p</span><span class="p">,</span> <span class="n">q</span> <span class="o">=</span> <span class="n">bernfrac</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">4</span><span class="p">)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">print</span><span class="p">(</span><span class="n">q</span><span class="p">)</span>
<span class="go">2338224387510</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">print</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="nb">str</span><span class="p">(</span><span class="n">p</span><span class="p">)))</span>
<span class="go">27692</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">print</span><span class="p">(</span><span class="n">mpf</span><span class="p">(</span><span class="n">p</span><span class="p">)</span> <span class="o">/</span> <span class="n">q</span><span class="p">)</span>
<span class="go">-9.04942396360948e+27677</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">print</span><span class="p">(</span><span class="n">bernoulli</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">4</span><span class="p">))</span>
<span class="go">-9.04942396360948e+27677</span>
</pre></div>
</div>
<div class="admonition note">
<p class="first admonition-title">Note</p>
<p class="last"><a title="mpmath.bernoulli" class="reference internal" href="#mpmath.bernoulli"><tt class="xref docutils literal"><span class="pre">bernoulli()</span></tt></a> computes a floating-point approximation
directly, without computing the exact fraction first.
This is much faster for large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>.</p>
</div>
<p><strong>Algorithm</strong></p>
<p><a title="mpmath.bernfrac" class="reference internal" href="#mpmath.bernfrac"><tt class="xref docutils literal"><span class="pre">bernfrac()</span></tt></a> works by computing the value of <img class="math" src="../_images/math/9be2b1de79fff66db67de3b18303bdd2912b7ae6.png" alt="B_n"/> numerically
and then using the von Staudt-Clausen theorem [1] to reconstruct
the exact fraction. For large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>, this is significantly faster than
computing <img class="math" src="../_images/math/1cd96ae0505c9d2838694d3a5b0aa6d11399258a.png" alt="B_1, B_2, \ldots, B_2"/> recursively with exact arithmetic.
The implementation has been tested for <img class="math" src="../_images/math/b68ef9ff01dd02bd64b149ecbb2066555a1df67e.png" alt="n = 10^m"/> up to <img class="math" src="../_images/math/970c67d43e2cd7f3e62f1e9afeb5d18508c8dd34.png" alt="m = 6"/>.</p>
<p>In practice, <a title="mpmath.bernfrac" class="reference internal" href="#mpmath.bernfrac"><tt class="xref docutils literal"><span class="pre">bernfrac()</span></tt></a> appears to be about three times
slower than the specialized program calcbn.exe [2]</p>
<p><strong>References</strong></p>
<ol class="arabic simple">
<li>MathWorld, von Staudt-Clausen Theorem:
<a class="reference external" href="http://mathworld.wolfram.com/vonStaudt-ClausenTheorem.html">http://mathworld.wolfram.com/vonStaudt-ClausenTheorem.html</a></li>
<li>The Bernoulli Number Page:
<a class="reference external" href="http://www.bernoulli.org/">http://www.bernoulli.org/</a></li>
</ol>
</dd></dl>

</div>
<div class="section" id="bernpoly">
<h3><tt class="xref docutils literal"><span class="pre">bernpoly()</span></tt><a class="headerlink" href="#bernpoly" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.bernpoly">
<tt class="descclassname">mpmath.</tt><tt class="descname">bernpoly</tt><big>(</big><em>n</em>, <em>z</em><big>)</big><a class="headerlink" href="#mpmath.bernpoly" title="Permalink to this definition">¶</a></dt>
<dd><p>Evaluates the Bernoulli polynomial <img class="math" src="../_images/math/d79269d871ce4723e28967c1f828f909535253d8.png" alt="B_n(z)"/>.</p>
<p>The first few Bernoulli polynomials are:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">6</span><span class="p">):</span>
<span class="gp">... </span>    <span class="n">nprint</span><span class="p">(</span><span class="n">chop</span><span class="p">(</span><span class="n">taylor</span><span class="p">(</span><span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">bernpoly</span><span class="p">(</span><span class="n">n</span><span class="p">,</span><span class="n">x</span><span class="p">),</span> <span class="mi">0</span><span class="p">,</span> <span class="n">n</span><span class="p">)))</span>
<span class="gp">...</span>
<span class="go">[1.0]</span>
<span class="go">[-0.5, 1.0]</span>
<span class="go">[0.166667, -1.0, 1.0]</span>
<span class="go">[0.0, 0.5, -1.5, 1.0]</span>
<span class="go">[-0.0333333, 0.0, 1.0, -2.0, 1.0]</span>
<span class="go">[0.0, -0.166667, 0.0, 1.66667, -2.5, 1.0]</span>
</pre></div>
</div>
<p>At <img class="math" src="../_images/math/aa9f0d97d2f39f78f05e05da40bf04f5a7c0726c.png" alt="z = 0"/>, the Bernoulli polynomial evaluates to a
Bernoulli number (see <a title="mpmath.bernoulli" class="reference internal" href="#mpmath.bernoulli"><tt class="xref docutils literal"><span class="pre">bernoulli()</span></tt></a>):</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">bernpoly</span><span class="p">(</span><span class="mi">12</span><span class="p">,</span> <span class="mi">0</span><span class="p">),</span> <span class="n">bernoulli</span><span class="p">(</span><span class="mi">12</span><span class="p">)</span>
<span class="go">(-0.253113553113553, -0.253113553113553)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bernpoly</span><span class="p">(</span><span class="mi">13</span><span class="p">,</span> <span class="mi">0</span><span class="p">),</span> <span class="n">bernoulli</span><span class="p">(</span><span class="mi">13</span><span class="p">)</span>
<span class="go">(0.0, 0.0)</span>
</pre></div>
</div>
<p>Evaluation is accurate for large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> and small <img class="math" src="../_images/math/b13f21416d84e13708696f34dea81026cda583c9.png" alt="z"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bernpoly</span><span class="p">(</span><span class="mi">100</span><span class="p">,</span> <span class="mf">0.5</span><span class="p">)</span>
<span class="go">2.838224957069370695926416e+78</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bernpoly</span><span class="p">(</span><span class="mi">1000</span><span class="p">,</span> <span class="mf">10.5</span><span class="p">)</span>
<span class="go">5.318704469415522036482914e+1769</span>
</pre></div>
</div>
</dd></dl>

</div>
</div>
<div class="section" id="euler-numbers-and-polynomials">
<h2>Euler numbers and polynomials<a class="headerlink" href="#euler-numbers-and-polynomials" title="Permalink to this headline">¶</a></h2>
<div class="section" id="eulernum">
<h3><tt class="xref docutils literal"><span class="pre">eulernum()</span></tt><a class="headerlink" href="#eulernum" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.eulernum">
<tt class="descclassname">mpmath.</tt><tt class="descname">eulernum</tt><big>(</big><em>n</em><big>)</big><a class="headerlink" href="#mpmath.eulernum" title="Permalink to this definition">¶</a></dt>
<dd><p>Gives the <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>-th Euler number, defined as the <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>-th derivative of
<img class="math" src="../_images/math/a54481a754624cd6f6aa6d88f130415905a83004.png" alt="\mathrm{sech}(t) = 1/\cosh(t)"/> evaluated at <img class="math" src="../_images/math/288e2a39b58d673182c57dc6214b702a448341ea.png" alt="t = 0"/>. Equivalently, the
Euler numbers give the coefficients of the Taylor series</p>
<div class="math">
<p><img src="../_images/math/3d79f8156c63633e6dfb072d38f3cf045a5381df.png" alt="\mathrm{sech}(t) = \sum_{n=0}^{\infty} \frac{E_n}{n!} t^n." /></p>
</div><p>The Euler numbers are closely related to Bernoulli numbers
and Bernoulli polynomials. They can also be evaluated in terms of
Euler polynomials (see <a title="mpmath.eulerpoly" class="reference internal" href="#mpmath.eulerpoly"><tt class="xref docutils literal"><span class="pre">eulerpoly()</span></tt></a>) as <img class="math" src="../_images/math/f9912f0eae7855e063b3980ce83501b953668bdf.png" alt="E_n = 2^n E_n(1/2)"/>.</p>
<p><strong>Examples</strong></p>
<p>Computing the first few Euler numbers and verifying that they
agree with the Taylor series:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="p">[</span><span class="n">eulernum</span><span class="p">(</span><span class="n">n</span><span class="p">)</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">11</span><span class="p">)]</span>
<span class="go">[1.0, 0.0, -1.0, 0.0, 5.0, 0.0, -61.0, 0.0, 1385.0, 0.0, -50521.0]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">chop</span><span class="p">(</span><span class="n">diffs</span><span class="p">(</span><span class="n">sech</span><span class="p">,</span> <span class="mi">0</span><span class="p">,</span> <span class="mi">10</span><span class="p">))</span>
<span class="go">[1.0, 0.0, -1.0, 0.0, 5.0, 0.0, -61.0, 0.0, 1385.0, 0.0, -50521.0]</span>
</pre></div>
</div>
<p>Euler numbers grow very rapidly. <a title="mpmath.eulernum" class="reference internal" href="#mpmath.eulernum"><tt class="xref docutils literal"><span class="pre">eulernum()</span></tt></a> efficiently
computes numerical approximations for large indices:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">eulernum</span><span class="p">(</span><span class="mi">50</span><span class="p">)</span>
<span class="go">-6.053285248188621896314384e+54</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulernum</span><span class="p">(</span><span class="mi">1000</span><span class="p">)</span>
<span class="go">3.887561841253070615257336e+2371</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulernum</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">20</span><span class="p">)</span>
<span class="go">4.346791453661149089338186e+1936958564106659551331</span>
</pre></div>
</div>
<p>Comparing with an asymptotic formula for the Euler numbers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">n</span> <span class="o">=</span> <span class="mi">10</span><span class="o">**</span><span class="mi">5</span>
<span class="gp">&gt;&gt;&gt; </span><span class="p">(</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span><span class="o">**</span><span class="p">(</span><span class="n">n</span><span class="o">//</span><span class="mi">2</span><span class="p">)</span> <span class="o">*</span> <span class="mi">8</span> <span class="o">*</span> <span class="n">sqrt</span><span class="p">(</span><span class="n">n</span><span class="o">/</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="p">))</span> <span class="o">*</span> <span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">n</span><span class="o">/</span><span class="p">(</span><span class="n">pi</span><span class="o">*</span><span class="n">e</span><span class="p">))</span><span class="o">**</span><span class="n">n</span>
<span class="go">3.69919063017432362805663e+436961</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulernum</span><span class="p">(</span><span class="n">n</span><span class="p">)</span>
<span class="go">3.699193712834466537941283e+436961</span>
</pre></div>
</div>
<p>Pass <tt class="docutils literal"><span class="pre">exact=True</span></tt> to obtain exact values of Euler numbers as integers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="k">print</span><span class="p">(</span><span class="n">eulernum</span><span class="p">(</span><span class="mi">50</span><span class="p">,</span> <span class="n">exact</span><span class="o">=</span><span class="bp">True</span><span class="p">))</span>
<span class="go">-6053285248188621896314383785111649088103498225146815121</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">print</span><span class="p">(</span><span class="n">eulernum</span><span class="p">(</span><span class="mi">200</span><span class="p">,</span> <span class="n">exact</span><span class="o">=</span><span class="bp">True</span><span class="p">)</span> <span class="o">%</span> <span class="mi">10</span><span class="o">**</span><span class="mi">10</span><span class="p">)</span>
<span class="go">1925859625</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulernum</span><span class="p">(</span><span class="mi">1001</span><span class="p">,</span> <span class="n">exact</span><span class="o">=</span><span class="bp">True</span><span class="p">)</span>
<span class="go">0</span>
</pre></div>
</div>
</dd></dl>

</div>
<div class="section" id="eulerpoly">
<h3><tt class="xref docutils literal"><span class="pre">eulerpoly()</span></tt><a class="headerlink" href="#eulerpoly" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.eulerpoly">
<tt class="descclassname">mpmath.</tt><tt class="descname">eulerpoly</tt><big>(</big><em>n</em>, <em>z</em><big>)</big><a class="headerlink" href="#mpmath.eulerpoly" title="Permalink to this definition">¶</a></dt>
<dd><p>Evaluates the Euler polynomial <img class="math" src="../_images/math/293ac7f9abdef50b197d8a3aee3acdbd1476db34.png" alt="E_n(z)"/>, defined by the generating function
representation</p>
<div class="math">
<p><img src="../_images/math/92cb5d84a9731bdcee624b4470f7e6f7002462c8.png" alt="\frac{2e^{zt}}{e^t+1} = \sum_{n=0}^\infty E_n(z) \frac{t^n}{n!}." /></p>
</div><p>The Euler polynomials may also be represented in terms of
Bernoulli polynomials (see <a title="mpmath.bernpoly" class="reference internal" href="#mpmath.bernpoly"><tt class="xref docutils literal"><span class="pre">bernpoly()</span></tt></a>) using various formulas, for
example</p>
<div class="math">
<p><img src="../_images/math/fc068cd19f0cb359f3e7fc72d246a210764646e0.png" alt="E_n(z) = \frac{2}{n+1} \left(
    B_n(z)-2^{n+1}B_n\left(\frac{z}{2}\right)
\right)." /></p>
</div><p>Special values include the Euler numbers <img class="math" src="../_images/math/f9912f0eae7855e063b3980ce83501b953668bdf.png" alt="E_n = 2^n E_n(1/2)"/> (see
<a title="mpmath.eulernum" class="reference internal" href="#mpmath.eulernum"><tt class="xref docutils literal"><span class="pre">eulernum()</span></tt></a>).</p>
<p><strong>Examples</strong></p>
<p>Computing the coefficients of the first few Euler polynomials:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">6</span><span class="p">):</span>
<span class="gp">... </span>    <span class="n">chop</span><span class="p">(</span><span class="n">taylor</span><span class="p">(</span><span class="k">lambda</span> <span class="n">z</span><span class="p">:</span> <span class="n">eulerpoly</span><span class="p">(</span><span class="n">n</span><span class="p">,</span><span class="n">z</span><span class="p">),</span> <span class="mi">0</span><span class="p">,</span> <span class="n">n</span><span class="p">))</span>
<span class="gp">...</span>
<span class="go">[1.0]</span>
<span class="go">[-0.5, 1.0]</span>
<span class="go">[0.0, -1.0, 1.0]</span>
<span class="go">[0.25, 0.0, -1.5, 1.0]</span>
<span class="go">[0.0, 1.0, 0.0, -2.0, 1.0]</span>
<span class="go">[-0.5, 0.0, 2.5, 0.0, -2.5, 1.0]</span>
</pre></div>
</div>
<p>Evaluation for arbitrary <img class="math" src="../_images/math/b13f21416d84e13708696f34dea81026cda583c9.png" alt="z"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">2</span><span class="p">,</span><span class="mi">3</span><span class="p">)</span>
<span class="go">6.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">5</span><span class="p">,</span><span class="mi">4</span><span class="p">)</span>
<span class="go">423.5</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">35</span><span class="p">,</span> <span class="mi">11111111112</span><span class="p">)</span>
<span class="go">3.994957561486776072734601e+351</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">4</span><span class="p">,</span> <span class="mi">10</span><span class="o">+</span><span class="mi">20</span><span class="n">j</span><span class="p">)</span>
<span class="go">(-47990.0 - 235980.0j)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">2</span><span class="p">,</span> <span class="s">&#39;-3.5e-5&#39;</span><span class="p">)</span>
<span class="go">0.000035001225</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">3</span><span class="p">,</span> <span class="mf">0.5</span><span class="p">)</span>
<span class="go">0.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">55</span><span class="p">,</span> <span class="o">-</span><span class="mi">10</span><span class="o">**</span><span class="mi">80</span><span class="p">)</span>
<span class="go">-1.0e+4400</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">5</span><span class="p">,</span> <span class="o">-</span><span class="n">inf</span><span class="p">)</span>
<span class="go">-inf</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span> <span class="o">-</span><span class="n">inf</span><span class="p">)</span>
<span class="go">+inf</span>
</pre></div>
</div>
<p>Computing Euler numbers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="mi">2</span><span class="o">**</span><span class="mi">26</span> <span class="o">*</span> <span class="n">eulerpoly</span><span class="p">(</span><span class="mi">26</span><span class="p">,</span><span class="mf">0.5</span><span class="p">)</span>
<span class="go">-4087072509293123892361.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulernum</span><span class="p">(</span><span class="mi">26</span><span class="p">)</span>
<span class="go">-4087072509293123892361.0</span>
</pre></div>
</div>
<p>Evaluation is accurate for large <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> and small <img class="math" src="../_images/math/b13f21416d84e13708696f34dea81026cda583c9.png" alt="z"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">100</span><span class="p">,</span> <span class="mf">0.5</span><span class="p">)</span>
<span class="go">2.29047999988194114177943e+108</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">1000</span><span class="p">,</span> <span class="mf">10.5</span><span class="p">)</span>
<span class="go">3.628120031122876847764566e+2070</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">eulerpoly</span><span class="p">(</span><span class="mi">10000</span><span class="p">,</span> <span class="mf">10.5</span><span class="p">)</span>
<span class="go">1.149364285543783412210773e+30688</span>
</pre></div>
</div>
</dd></dl>

</div>
</div>
<div class="section" id="bell-numbers-and-polynomials">
<h2>Bell numbers and polynomials<a class="headerlink" href="#bell-numbers-and-polynomials" title="Permalink to this headline">¶</a></h2>
<div class="section" id="bell">
<h3><tt class="xref docutils literal"><span class="pre">bell()</span></tt><a class="headerlink" href="#bell" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.bell">
<tt class="descclassname">mpmath.</tt><tt class="descname">bell</tt><big>(</big><em>n</em>, <em>x</em><big>)</big><a class="headerlink" href="#mpmath.bell" title="Permalink to this definition">¶</a></dt>
<dd><p>For <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> a nonnegative integer, <tt class="docutils literal"><span class="pre">bell(n,x)</span></tt> evaluates the Bell
polynomial <img class="math" src="../_images/math/6343a0dfcc1bb7d992dcdf886dc7c3947e1a3239.png" alt="B_n(x)"/>, the first few of which are</p>
<div class="math">
<p><img src="../_images/math/370eaaa999927a40127739b601e85935724735d3.png" alt="B_0(x) = 1

B_1(x) = x

B_2(x) = x^2+x

B_3(x) = x^3+3x^2+x" /></p>
</div><p>If <img class="math" src="../_images/math/34310f2f36ed6d724838edb08788ee62acb33386.png" alt="x = 1"/> or <a title="mpmath.bell" class="reference internal" href="#mpmath.bell"><tt class="xref docutils literal"><span class="pre">bell()</span></tt></a> is called with only one argument, it
gives the <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>-th Bell number <img class="math" src="../_images/math/9be2b1de79fff66db67de3b18303bdd2912b7ae6.png" alt="B_n"/>, which is the number of
partitions of a set with <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/> elements. By setting the precision to
at least <img class="math" src="../_images/math/816b0b85300fda85bdb84aa1858eadc31859deed.png" alt="\log_{10} B_n"/> digits, <a title="mpmath.bell" class="reference internal" href="#mpmath.bell"><tt class="xref docutils literal"><span class="pre">bell()</span></tt></a> provides fast
calculation of exact Bell numbers.</p>
<p>In general, <a title="mpmath.bell" class="reference internal" href="#mpmath.bell"><tt class="xref docutils literal"><span class="pre">bell()</span></tt></a> computes</p>
<div class="math">
<p><img src="../_images/math/14f27b46bfc53a41da9762cad023072e45efb5a0.png" alt="B_n(x) = e^{-x} \left(\mathrm{sinc}(\pi n) + E_n(x)\right)" /></p>
</div><p>where <img class="math" src="../_images/math/65e3e61fce7353dfe4ed8461c9c994afc38c2492.png" alt="E_n(x)"/> is the generalized exponential function implemented
by <a title="mpmath.polyexp" class="reference external" href="zeta.html#mpmath.polyexp"><tt class="xref docutils literal"><span class="pre">polyexp()</span></tt></a>. This is an extension of Dobinski&#8217;s formula [1],
where the modification is the sinc term ensuring that <img class="math" src="../_images/math/6343a0dfcc1bb7d992dcdf886dc7c3947e1a3239.png" alt="B_n(x)"/> is
continuous in <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>; <a title="mpmath.bell" class="reference internal" href="#mpmath.bell"><tt class="xref docutils literal"><span class="pre">bell()</span></tt></a> can thus be evaluated,
differentiated, etc for arbitrary complex arguments.</p>
<p><strong>Examples</strong></p>
<p>Simple evaluations:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">0</span><span class="p">,</span> <span class="mf">2.5</span><span class="p">)</span>
<span class="go">1.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">1</span><span class="p">,</span> <span class="mf">2.5</span><span class="p">)</span>
<span class="go">2.5</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">2</span><span class="p">,</span> <span class="mf">2.5</span><span class="p">)</span>
<span class="go">8.75</span>
</pre></div>
</div>
<p>Evaluation for arbitrary complex arguments:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mf">5.75</span><span class="o">+</span><span class="mi">1</span><span class="n">j</span><span class="p">,</span> <span class="mi">2</span><span class="o">-</span><span class="mi">3</span><span class="n">j</span><span class="p">)</span>
<span class="go">(-10767.71345136587098445143 - 15449.55065599872579097221j)</span>
</pre></div>
</div>
<p>The first few Bell polynomials:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">k</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">7</span><span class="p">):</span>
<span class="gp">... </span>    <span class="n">nprint</span><span class="p">(</span><span class="n">taylor</span><span class="p">(</span><span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">bell</span><span class="p">(</span><span class="n">k</span><span class="p">,</span><span class="n">x</span><span class="p">),</span> <span class="mi">0</span><span class="p">,</span> <span class="n">k</span><span class="p">))</span>
<span class="gp">...</span>
<span class="go">[1.0]</span>
<span class="go">[0.0, 1.0]</span>
<span class="go">[0.0, 1.0, 1.0]</span>
<span class="go">[0.0, 1.0, 3.0, 1.0]</span>
<span class="go">[0.0, 1.0, 7.0, 6.0, 1.0]</span>
<span class="go">[0.0, 1.0, 15.0, 25.0, 10.0, 1.0]</span>
<span class="go">[0.0, 1.0, 31.0, 90.0, 65.0, 15.0, 1.0]</span>
</pre></div>
</div>
<p>The first few Bell numbers and complementary Bell numbers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="p">[</span><span class="nb">int</span><span class="p">(</span><span class="n">bell</span><span class="p">(</span><span class="n">k</span><span class="p">))</span> <span class="k">for</span> <span class="n">k</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">10</span><span class="p">)]</span>
<span class="go">[1, 1, 2, 5, 15, 52, 203, 877, 4140, 21147]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="p">[</span><span class="nb">int</span><span class="p">(</span><span class="n">bell</span><span class="p">(</span><span class="n">k</span><span class="p">,</span><span class="o">-</span><span class="mi">1</span><span class="p">))</span> <span class="k">for</span> <span class="n">k</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">10</span><span class="p">)]</span>
<span class="go">[1, -1, 0, 1, 1, -2, -9, -9, 50, 267]</span>
</pre></div>
</div>
<p>Large Bell numbers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">50</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">50</span><span class="p">)</span>
<span class="go">185724268771078270438257767181908917499221852770.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">50</span><span class="p">,</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span>
<span class="go">-29113173035759403920216141265491160286912.0</span>
</pre></div>
</div>
<p>Some even larger values:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">1000</span><span class="p">,</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span>
<span class="go">-1.237132026969293954162816e+1869</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">1000</span><span class="p">)</span>
<span class="go">2.989901335682408421480422e+1927</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">1000</span><span class="p">,</span><span class="mi">2</span><span class="p">)</span>
<span class="go">6.591553486811969380442171e+1987</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">bell</span><span class="p">(</span><span class="mi">1000</span><span class="p">,</span><span class="mf">100.5</span><span class="p">)</span>
<span class="go">9.101014101401543575679639e+2529</span>
</pre></div>
</div>
<p>A determinant identity satisfied by Bell numbers:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">N</span> <span class="o">=</span> <span class="mi">8</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">det</span><span class="p">([[</span><span class="n">bell</span><span class="p">(</span><span class="n">k</span><span class="o">+</span><span class="n">j</span><span class="p">)</span> <span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="n">N</span><span class="p">)]</span> <span class="k">for</span> <span class="n">k</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="n">N</span><span class="p">)])</span>
<span class="go">125411328000.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">superfac</span><span class="p">(</span><span class="n">N</span><span class="o">-</span><span class="mi">1</span><span class="p">)</span>
<span class="go">125411328000.0</span>
</pre></div>
</div>
<p><strong>References</strong></p>
<ol class="arabic simple">
<li><a class="reference external" href="http://mathworld.wolfram.com/DobinskisFormula.html">http://mathworld.wolfram.com/DobinskisFormula.html</a></li>
</ol>
</dd></dl>

</div>
</div>
<div class="section" id="prime-counting-functions">
<h2>Prime counting functions<a class="headerlink" href="#prime-counting-functions" title="Permalink to this headline">¶</a></h2>
<div class="section" id="primepi">
<h3><tt class="xref docutils literal"><span class="pre">primepi()</span></tt><a class="headerlink" href="#primepi" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.primepi">
<tt class="descclassname">mpmath.</tt><tt class="descname">primepi</tt><big>(</big><em>x</em><big>)</big><a class="headerlink" href="#mpmath.primepi" title="Permalink to this definition">¶</a></dt>
<dd><p>Evaluates the prime counting function, <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/>, which gives
the number of primes less than or equal to <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/>. The argument
<img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/> may be fractional.</p>
<p>The prime counting function is very expensive to evaluate
precisely for large <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/>, and the present implementation is
not optimized in any way. For numerical approximation of the
prime counting function, it is better to use <a title="mpmath.primepi2" class="reference internal" href="#mpmath.primepi2"><tt class="xref docutils literal"><span class="pre">primepi2()</span></tt></a>
or <a title="mpmath.riemannr" class="reference internal" href="#mpmath.riemannr"><tt class="xref docutils literal"><span class="pre">riemannr()</span></tt></a>.</p>
<p>Some values of the prime counting function:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="p">[</span><span class="n">primepi</span><span class="p">(</span><span class="n">k</span><span class="p">)</span> <span class="k">for</span> <span class="n">k</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">20</span><span class="p">)]</span>
<span class="go">[0, 0, 1, 2, 2, 3, 3, 4, 4, 4, 4, 5, 5, 6, 6, 6, 6, 7, 7, 8]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi</span><span class="p">(</span><span class="mf">3.5</span><span class="p">)</span>
<span class="go">2</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi</span><span class="p">(</span><span class="mi">100000</span><span class="p">)</span>
<span class="go">9592</span>
</pre></div>
</div>
</dd></dl>

</div>
<div class="section" id="primepi2">
<h3><tt class="xref docutils literal"><span class="pre">primepi2()</span></tt><a class="headerlink" href="#primepi2" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.primepi2">
<tt class="descclassname">mpmath.</tt><tt class="descname">primepi2</tt><big>(</big><em>x</em><big>)</big><a class="headerlink" href="#mpmath.primepi2" title="Permalink to this definition">¶</a></dt>
<dd><p>Returns an interval (as an <tt class="docutils literal"><span class="pre">mpi</span></tt> instance) providing bounds
for the value of the prime counting function <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/>. For small
<img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/>, <a title="mpmath.primepi2" class="reference internal" href="#mpmath.primepi2"><tt class="xref docutils literal"><span class="pre">primepi2()</span></tt></a> returns an exact interval based on
the output of <a title="mpmath.primepi" class="reference internal" href="#mpmath.primepi"><tt class="xref docutils literal"><span class="pre">primepi()</span></tt></a>. For <img class="math" src="../_images/math/f54fa98073b37791da5e80d6e21a5282536adf29.png" alt="x &gt; 2656"/>, a loose interval
based on Schoenfeld&#8217;s inequality</p>
<div class="math">
<p><img src="../_images/math/ec56afd657455289103a265e4b38ae8ec1a635b4.png" alt="|\pi(x) - \mathrm{li}(x)| &lt; \frac{\sqrt x \log x}{8 \pi}" /></p>
</div><p>is returned. This estimate is rigorous assuming the truth of
the Riemann hypothesis, and can be computed very quickly.</p>
<p><strong>Examples</strong></p>
<p>Exact values of the prime counting function for small <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">iv</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">iv</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi2</span><span class="p">(</span><span class="mi">10</span><span class="p">)</span>
<span class="go">[4.0, 4.0]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi2</span><span class="p">(</span><span class="mi">100</span><span class="p">)</span>
<span class="go">[25.0, 25.0]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi2</span><span class="p">(</span><span class="mi">1000</span><span class="p">)</span>
<span class="go">[168.0, 168.0]</span>
</pre></div>
</div>
<p>Loose intervals are generated for moderately large <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">primepi2</span><span class="p">(</span><span class="mi">10000</span><span class="p">),</span> <span class="n">primepi</span><span class="p">(</span><span class="mi">10000</span><span class="p">)</span>
<span class="go">([1209.0, 1283.0], 1229)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi2</span><span class="p">(</span><span class="mi">50000</span><span class="p">),</span> <span class="n">primepi</span><span class="p">(</span><span class="mi">50000</span><span class="p">)</span>
<span class="go">([5070.0, 5263.0], 5133)</span>
</pre></div>
</div>
<p>As <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/> increases, the absolute error gets worse while the relative
error improves. The exact value of <img class="math" src="../_images/math/11c7d49283e43405ab297a52f1f2aa0de7d32894.png" alt="\pi(10^{23})"/> is
1925320391606803968923, and <a title="mpmath.primepi2" class="reference internal" href="#mpmath.primepi2"><tt class="xref docutils literal"><span class="pre">primepi2()</span></tt></a> gives 9 significant
digits:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">p</span> <span class="o">=</span> <span class="n">primepi2</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">23</span><span class="p">)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">p</span>
<span class="go">[1.9253203909477020467e+21, 1.925320392280406229e+21]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mpf</span><span class="p">(</span><span class="n">p</span><span class="o">.</span><span class="n">delta</span><span class="p">)</span> <span class="o">/</span> <span class="n">mpf</span><span class="p">(</span><span class="n">p</span><span class="o">.</span><span class="n">a</span><span class="p">)</span>
<span class="go">6.9219865355293e-10</span>
</pre></div>
</div>
<p>A more precise, nonrigorous estimate for <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/> can be
obtained using the Riemann R function (<a title="mpmath.riemannr" class="reference internal" href="#mpmath.riemannr"><tt class="xref docutils literal"><span class="pre">riemannr()</span></tt></a>).
For large enough <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/>, the value returned by <a title="mpmath.primepi2" class="reference internal" href="#mpmath.primepi2"><tt class="xref docutils literal"><span class="pre">primepi2()</span></tt></a>
essentially amounts to a small perturbation of the value returned by
<a title="mpmath.riemannr" class="reference internal" href="#mpmath.riemannr"><tt class="xref docutils literal"><span class="pre">riemannr()</span></tt></a>:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">primepi2</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">100</span><span class="p">)</span>
<span class="go">[4.3619719871407024816e+97, 4.3619719871407032404e+97]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">riemannr</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">100</span><span class="p">)</span>
<span class="go">4.3619719871407e+97</span>
</pre></div>
</div>
</dd></dl>

</div>
<div class="section" id="riemannr">
<h3><tt class="xref docutils literal"><span class="pre">riemannr()</span></tt><a class="headerlink" href="#riemannr" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.riemannr">
<tt class="descclassname">mpmath.</tt><tt class="descname">riemannr</tt><big>(</big><em>x</em><big>)</big><a class="headerlink" href="#mpmath.riemannr" title="Permalink to this definition">¶</a></dt>
<dd><p>Evaluates the Riemann R function, a smooth approximation of the
prime counting function <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/> (see <a title="mpmath.primepi" class="reference internal" href="#mpmath.primepi"><tt class="xref docutils literal"><span class="pre">primepi()</span></tt></a>). The Riemann
R function gives a fast numerical approximation useful e.g. to
roughly estimate the number of primes in a given interval.</p>
<p>The Riemann R function is computed using the rapidly convergent Gram
series,</p>
<div class="math">
<p><img src="../_images/math/9abdb37e590968789986e9be89fd6487c1dbd72b.png" alt="R(x) = 1 + \sum_{k=1}^{\infty}
    \frac{\log^k x}{k k! \zeta(k+1)}." /></p>
</div><p>From the Gram series, one sees that the Riemann R function is a
well-defined analytic function (except for a branch cut along
the negative real half-axis); it can be evaluated for arbitrary
real or complex arguments.</p>
<p>The Riemann R function gives a very accurate approximation
of the prime counting function. For example, it is wrong by at
most 2 for <img class="math" src="../_images/math/9c7c4e49985884250ef23345752e698d917257d1.png" alt="x &lt; 1000"/>, and for <img class="math" src="../_images/math/5135063eaa0d54a454c2af8e685b4af83222a871.png" alt="x = 10^9"/> differs from the exact
value of <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/> by 79, or less than two parts in a million.
It is about 10 times more accurate than the logarithmic integral
estimate (see <a title="mpmath.li" class="reference external" href="expintegrals.html#mpmath.li"><tt class="xref docutils literal"><span class="pre">li()</span></tt></a>), which however is even faster to evaluate.
It is orders of magnitude more accurate than the extremely
fast <img class="math" src="../_images/math/a07bb8ab7d72be94874de0133e8353cd654b9a11.png" alt="x/\log x"/> estimate.</p>
<p><strong>Examples</strong></p>
<p>For small arguments, the Riemann R function almost exactly
gives the prime counting function if rounded to the nearest
integer:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">primepi</span><span class="p">(</span><span class="mi">50</span><span class="p">),</span> <span class="n">riemannr</span><span class="p">(</span><span class="mi">50</span><span class="p">)</span>
<span class="go">(15, 14.9757023241462)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="nb">max</span><span class="p">(</span><span class="nb">abs</span><span class="p">(</span><span class="n">primepi</span><span class="p">(</span><span class="n">n</span><span class="p">)</span><span class="o">-</span><span class="nb">int</span><span class="p">(</span><span class="nb">round</span><span class="p">(</span><span class="n">riemannr</span><span class="p">(</span><span class="n">n</span><span class="p">))))</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">100</span><span class="p">))</span>
<span class="go">1</span>
<span class="gp">&gt;&gt;&gt; </span><span class="nb">max</span><span class="p">(</span><span class="nb">abs</span><span class="p">(</span><span class="n">primepi</span><span class="p">(</span><span class="n">n</span><span class="p">)</span><span class="o">-</span><span class="nb">int</span><span class="p">(</span><span class="nb">round</span><span class="p">(</span><span class="n">riemannr</span><span class="p">(</span><span class="n">n</span><span class="p">))))</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">300</span><span class="p">))</span>
<span class="go">2</span>
</pre></div>
</div>
<p>The Riemann R function can be evaluated for arguments far too large
for exact determination of <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/> to be computationally
feasible with any presently known algorithm:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">riemannr</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">30</span><span class="p">)</span>
<span class="go">1.46923988977204e+28</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">riemannr</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">100</span><span class="p">)</span>
<span class="go">4.3619719871407e+97</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">riemannr</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="mi">1000</span><span class="p">)</span>
<span class="go">4.3448325764012e+996</span>
</pre></div>
</div>
<p>A comparison of the Riemann R function and logarithmic integral estimates
for <img class="math" src="../_images/math/bf5eda8e931643635bbb49a1b8780a45fcba891e.png" alt="\pi(x)"/> using exact values of <img class="math" src="../_images/math/74c776a0c3180f43dd1ce02b63c33af5a9568dc7.png" alt="\pi(10^n)"/> up to <img class="math" src="../_images/math/1819f547a53f667e1a44a49ba2bec08c28d7f222.png" alt="n = 9"/>.
The fractional error is shown in parentheses:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">exact</span> <span class="o">=</span> <span class="p">[</span><span class="mi">4</span><span class="p">,</span><span class="mi">25</span><span class="p">,</span><span class="mi">168</span><span class="p">,</span><span class="mi">1229</span><span class="p">,</span><span class="mi">9592</span><span class="p">,</span><span class="mi">78498</span><span class="p">,</span><span class="mi">664579</span><span class="p">,</span><span class="mi">5761455</span><span class="p">,</span><span class="mi">50847534</span><span class="p">]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">n</span><span class="p">,</span> <span class="n">p</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">exact</span><span class="p">):</span>
<span class="gp">... </span>    <span class="n">n</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="gp">... </span>    <span class="n">r</span><span class="p">,</span> <span class="n">l</span> <span class="o">=</span> <span class="n">riemannr</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="n">n</span><span class="p">),</span> <span class="n">li</span><span class="p">(</span><span class="mi">10</span><span class="o">**</span><span class="n">n</span><span class="p">)</span>
<span class="gp">... </span>    <span class="n">rerr</span><span class="p">,</span> <span class="n">lerr</span> <span class="o">=</span> <span class="n">nstr</span><span class="p">((</span><span class="n">r</span><span class="o">-</span><span class="n">p</span><span class="p">)</span><span class="o">/</span><span class="n">p</span><span class="p">,</span><span class="mi">3</span><span class="p">),</span> <span class="n">nstr</span><span class="p">((</span><span class="n">l</span><span class="o">-</span><span class="n">p</span><span class="p">)</span><span class="o">/</span><span class="n">p</span><span class="p">,</span><span class="mi">3</span><span class="p">)</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="s">&quot;</span><span class="si">%i</span><span class="s"> </span><span class="si">%i</span><span class="s"> </span><span class="si">%s</span><span class="s">(</span><span class="si">%s</span><span class="s">) </span><span class="si">%s</span><span class="s">(</span><span class="si">%s</span><span class="s">)&quot;</span> <span class="o">%</span> <span class="p">(</span><span class="n">n</span><span class="p">,</span> <span class="n">p</span><span class="p">,</span> <span class="n">r</span><span class="p">,</span> <span class="n">rerr</span><span class="p">,</span> <span class="n">l</span><span class="p">,</span> <span class="n">lerr</span><span class="p">))</span>
<span class="gp">...</span>
<span class="go">1 4 4.56458314100509(0.141) 6.1655995047873(0.541)</span>
<span class="go">2 25 25.6616332669242(0.0265) 30.1261415840796(0.205)</span>
<span class="go">3 168 168.359446281167(0.00214) 177.609657990152(0.0572)</span>
<span class="go">4 1229 1226.93121834343(-0.00168) 1246.13721589939(0.0139)</span>
<span class="go">5 9592 9587.43173884197(-0.000476) 9629.8090010508(0.00394)</span>
<span class="go">6 78498 78527.3994291277(0.000375) 78627.5491594622(0.00165)</span>
<span class="go">7 664579 664667.447564748(0.000133) 664918.405048569(0.000511)</span>
<span class="go">8 5761455 5761551.86732017(1.68e-5) 5762209.37544803(0.000131)</span>
<span class="go">9 50847534 50847455.4277214(-1.55e-6) 50849234.9570018(3.35e-5)</span>
</pre></div>
</div>
<p>The derivative of the Riemann R function gives the approximate
probability for a number of magnitude <img class="math" src="../_images/math/26eeb5258ca5099acf8fe96b2a1049c48c89a5e6.png" alt="x"/> to be prime:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">diff</span><span class="p">(</span><span class="n">riemannr</span><span class="p">,</span> <span class="mi">1000</span><span class="p">)</span>
<span class="go">0.141903028110784</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mpf</span><span class="p">(</span><span class="n">primepi</span><span class="p">(</span><span class="mi">1050</span><span class="p">)</span> <span class="o">-</span> <span class="n">primepi</span><span class="p">(</span><span class="mi">950</span><span class="p">))</span> <span class="o">/</span> <span class="mi">100</span>
<span class="go">0.15</span>
</pre></div>
</div>
<p>Evaluation is supported for arbitrary arguments and at arbitrary
precision:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">30</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">riemannr</span><span class="p">(</span><span class="mf">7.5</span><span class="p">)</span>
<span class="go">3.72934743264966261918857135136</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">riemannr</span><span class="p">(</span><span class="o">-</span><span class="mi">4</span><span class="o">+</span><span class="mi">2</span><span class="n">j</span><span class="p">)</span>
<span class="go">(-0.551002208155486427591793957644 + 2.16966398138119450043195899746j)</span>
</pre></div>
</div>
</dd></dl>

</div>
</div>
<div class="section" id="cyclotomic-polynomials">
<h2>Cyclotomic polynomials<a class="headerlink" href="#cyclotomic-polynomials" title="Permalink to this headline">¶</a></h2>
<div class="section" id="cyclotomic">
<h3><tt class="xref docutils literal"><span class="pre">cyclotomic()</span></tt><a class="headerlink" href="#cyclotomic" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.cyclotomic">
<tt class="descclassname">mpmath.</tt><tt class="descname">cyclotomic</tt><big>(</big><em>n</em>, <em>x</em><big>)</big><a class="headerlink" href="#mpmath.cyclotomic" title="Permalink to this definition">¶</a></dt>
<dd><p>Evaluates the cyclotomic polynomial <img class="math" src="../_images/math/448d49f4d827bb1e6c7927b2235a00d1abebe582.png" alt="\Phi_n(x)"/>, defined by</p>
<div class="math">
<p><img src="../_images/math/d82a3061f68e10836bed951249787f439dcb1cdb.png" alt="\Phi_n(x) = \prod_{\zeta} (x - \zeta)" /></p>
</div><p>where <img class="math" src="../_images/math/fda38f45a64bbd4451aedeeb4584e01d91239172.png" alt="\zeta"/> ranges over all primitive <img class="math" src="../_images/math/174fadd07fd54c9afe288e96558c92e0c1da733a.png" alt="n"/>-th roots of unity
(see <a title="mpmath.unitroots" class="reference external" href="powers.html#mpmath.unitroots"><tt class="xref docutils literal"><span class="pre">unitroots()</span></tt></a>). An equivalent representation, used
for computation, is</p>
<div class="math">
<p><img src="../_images/math/da61a04df78eb35d4636f148b159cc74321b5879.png" alt="\Phi_n(x) = \prod_{d\mid n}(x^d-1)^{\mu(n/d)} = \Phi_n(x)" /></p>
</div><p>where <img class="math" src="../_images/math/ad8136a3946faf34a8e134e613b2b4089e91bfa9.png" alt="\mu(m)"/> denotes the Moebius function. The cyclotomic
polynomials are integer polynomials, the first of which can be
written explicitly as</p>
<div class="math">
<p><img src="../_images/math/41183cd1b6977cae3d96f939d5aebfa80c97c53f.png" alt="\Phi_0(x) = 1

\Phi_1(x) = x - 1

\Phi_2(x) = x + 1

\Phi_3(x) = x^3 + x^2 + 1

\Phi_4(x) = x^2 + 1

\Phi_5(x) = x^4 + x^3 + x^2 + x + 1

\Phi_6(x) = x^2 - x + 1" /></p>
</div><p><strong>Examples</strong></p>
<p>The coefficients of low-order cyclotomic polynomials can be recovered
using Taylor expansion:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">15</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">9</span><span class="p">):</span>
<span class="gp">... </span>    <span class="n">p</span> <span class="o">=</span> <span class="n">chop</span><span class="p">(</span><span class="n">taylor</span><span class="p">(</span><span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">cyclotomic</span><span class="p">(</span><span class="n">n</span><span class="p">,</span><span class="n">x</span><span class="p">),</span> <span class="mi">0</span><span class="p">,</span> <span class="mi">10</span><span class="p">))</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="s">&quot;</span><span class="si">%s</span><span class="s"> </span><span class="si">%s</span><span class="s">&quot;</span> <span class="o">%</span> <span class="p">(</span><span class="n">n</span><span class="p">,</span> <span class="n">nstr</span><span class="p">(</span><span class="n">p</span><span class="p">[:</span><span class="mi">10</span><span class="o">+</span><span class="mi">1</span><span class="o">-</span><span class="n">p</span><span class="p">[::</span><span class="o">-</span><span class="mi">1</span><span class="p">]</span><span class="o">.</span><span class="n">index</span><span class="p">(</span><span class="mi">1</span><span class="p">)])))</span>
<span class="gp">...</span>
<span class="go">0 [1.0]</span>
<span class="go">1 [-1.0, 1.0]</span>
<span class="go">2 [1.0, 1.0]</span>
<span class="go">3 [1.0, 1.0, 1.0]</span>
<span class="go">4 [1.0, 0.0, 1.0]</span>
<span class="go">5 [1.0, 1.0, 1.0, 1.0, 1.0]</span>
<span class="go">6 [1.0, -1.0, 1.0]</span>
<span class="go">7 [1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0]</span>
<span class="go">8 [1.0, 0.0, 0.0, 0.0, 1.0]</span>
</pre></div>
</div>
<p>The definition as a product over primitive roots may be checked
by computing the product explicitly (for a real argument, this
method will generally introduce numerical noise in the imaginary
part):</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">z</span> <span class="o">=</span> <span class="mi">3</span><span class="o">+</span><span class="mi">4</span><span class="n">j</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">cyclotomic</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span> <span class="n">z</span><span class="p">)</span>
<span class="go">(-419.0 - 360.0j)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fprod</span><span class="p">(</span><span class="n">z</span><span class="o">-</span><span class="n">r</span> <span class="k">for</span> <span class="n">r</span> <span class="ow">in</span> <span class="n">unitroots</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span> <span class="n">primitive</span><span class="o">=</span><span class="bp">True</span><span class="p">))</span>
<span class="go">(-419.0 - 360.0j)</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">z</span> <span class="o">=</span> <span class="mi">3</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">cyclotomic</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span> <span class="n">z</span><span class="p">)</span>
<span class="go">61.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fprod</span><span class="p">(</span><span class="n">z</span><span class="o">-</span><span class="n">r</span> <span class="k">for</span> <span class="n">r</span> <span class="ow">in</span> <span class="n">unitroots</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span> <span class="n">primitive</span><span class="o">=</span><span class="bp">True</span><span class="p">))</span>
<span class="go">(61.0 - 3.146045605088568607055454e-25j)</span>
</pre></div>
</div>
<p>Up to permutation, the roots of a given cyclotomic polynomial
can be checked to agree with the list of primitive roots:</p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="n">p</span> <span class="o">=</span> <span class="n">taylor</span><span class="p">(</span><span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">cyclotomic</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span><span class="n">x</span><span class="p">),</span> <span class="mi">0</span><span class="p">,</span> <span class="mi">6</span><span class="p">)[:</span><span class="mi">3</span><span class="p">]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">r</span> <span class="ow">in</span> <span class="n">polyroots</span><span class="p">(</span><span class="n">p</span><span class="p">[::</span><span class="o">-</span><span class="mi">1</span><span class="p">]):</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="n">r</span><span class="p">)</span>
<span class="gp">...</span>
<span class="go">(0.5 - 0.8660254037844386467637232j)</span>
<span class="go">(0.5 + 0.8660254037844386467637232j)</span>
<span class="go">&gt;&gt;&gt;</span>
<span class="gp">&gt;&gt;&gt; </span><span class="k">for</span> <span class="n">r</span> <span class="ow">in</span> <span class="n">unitroots</span><span class="p">(</span><span class="mi">6</span><span class="p">,</span> <span class="n">primitive</span><span class="o">=</span><span class="bp">True</span><span class="p">):</span>
<span class="gp">... </span>    <span class="k">print</span><span class="p">(</span><span class="n">r</span><span class="p">)</span>
<span class="gp">...</span>
<span class="go">(0.5 + 0.8660254037844386467637232j)</span>
<span class="go">(0.5 - 0.8660254037844386467637232j)</span>
</pre></div>
</div>
</dd></dl>

</div>
</div>
<div class="section" id="arithmetic-functions">
<h2>Arithmetic functions<a class="headerlink" href="#arithmetic-functions" title="Permalink to this headline">¶</a></h2>
<div class="section" id="mangoldt">
<h3><tt class="xref docutils literal"><span class="pre">mangoldt()</span></tt><a class="headerlink" href="#mangoldt" title="Permalink to this headline">¶</a></h3>
<dl class="function">
<dt id="mpmath.mangoldt">
<tt class="descclassname">mpmath.</tt><tt class="descname">mangoldt</tt><big>(</big><em>n</em><big>)</big><a class="headerlink" href="#mpmath.mangoldt" title="Permalink to this definition">¶</a></dt>
<dd><p>Evaluates the von Mangoldt function <img class="math" src="../_images/math/962d7bd72e1c61e68bafaacfcfdefc37e46e884d.png" alt="\Lambda(n) = \log p"/>
if <img class="math" src="../_images/math/4879efbc807f21af3d3a83d8a3c1dbac4800e229.png" alt="n = p^k"/> a power of a prime, and <img class="math" src="../_images/math/6174a5cd45901b03273c42033ed99cf936a516c7.png" alt="\Lambda(n) = 0"/> otherwise.</p>
<p><strong>Examples</strong></p>
<div class="highlight-python"><div class="highlight"><pre><span class="gp">&gt;&gt;&gt; </span><span class="kn">from</span> <span class="nn">mpmath</span> <span class="kn">import</span> <span class="o">*</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mp</span><span class="o">.</span><span class="n">dps</span> <span class="o">=</span> <span class="mi">25</span><span class="p">;</span> <span class="n">mp</span><span class="o">.</span><span class="n">pretty</span> <span class="o">=</span> <span class="bp">True</span>
<span class="gp">&gt;&gt;&gt; </span><span class="p">[</span><span class="n">mangoldt</span><span class="p">(</span><span class="n">n</span><span class="p">)</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="o">-</span><span class="mi">2</span><span class="p">,</span><span class="mi">3</span><span class="p">)]</span>
<span class="go">[0.0, 0.0, 0.0, 0.0, 0.6931471805599453094172321]</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mangoldt</span><span class="p">(</span><span class="mi">6</span><span class="p">)</span>
<span class="go">0.0</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mangoldt</span><span class="p">(</span><span class="mi">7</span><span class="p">)</span>
<span class="go">1.945910149055313305105353</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">mangoldt</span><span class="p">(</span><span class="mi">8</span><span class="p">)</span>
<span class="go">0.6931471805599453094172321</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fsum</span><span class="p">(</span><span class="n">mangoldt</span><span class="p">(</span><span class="n">n</span><span class="p">)</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">101</span><span class="p">))</span>
<span class="go">94.04531122935739224600493</span>
<span class="gp">&gt;&gt;&gt; </span><span class="n">fsum</span><span class="p">(</span><span class="n">mangoldt</span><span class="p">(</span><span class="n">n</span><span class="p">)</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">10001</span><span class="p">))</span>
<span class="go">10013.39669326311478372032</span>
</pre></div>
</div>
</dd></dl>

</div>
</div>
</div>


          </div>
        </div>
      </div>
      <div class="sphinxsidebar">
        <div class="sphinxsidebarwrapper">
            <h3><a href="../index.html">Table Of Contents</a></h3>
            <ul>
<li><a class="reference external" href="#">Number-theoretical, combinatorial and integer functions</a><ul>
<li><a class="reference external" href="#fibonacci-numbers">Fibonacci numbers</a><ul>
<li><a class="reference external" href="#fibonacci-fib"><tt class="docutils literal"><span class="pre">fibonacci()</span></tt>/<tt class="docutils literal"><span class="pre">fib()</span></tt></a></li>
</ul>
</li>
<li><a class="reference external" href="#bernoulli-numbers-and-polynomials">Bernoulli numbers and polynomials</a><ul>
<li><a class="reference external" href="#bernoulli"><tt class="docutils literal"><span class="pre">bernoulli()</span></tt></a></li>
<li><a class="reference external" href="#bernfrac"><tt class="docutils literal"><span class="pre">bernfrac()</span></tt></a></li>
<li><a class="reference external" href="#bernpoly"><tt class="docutils literal"><span class="pre">bernpoly()</span></tt></a></li>
</ul>
</li>
<li><a class="reference external" href="#euler-numbers-and-polynomials">Euler numbers and polynomials</a><ul>
<li><a class="reference external" href="#eulernum"><tt class="docutils literal"><span class="pre">eulernum()</span></tt></a></li>
<li><a class="reference external" href="#eulerpoly"><tt class="docutils literal"><span class="pre">eulerpoly()</span></tt></a></li>
</ul>
</li>
<li><a class="reference external" href="#bell-numbers-and-polynomials">Bell numbers and polynomials</a><ul>
<li><a class="reference external" href="#bell"><tt class="docutils literal"><span class="pre">bell()</span></tt></a></li>
</ul>
</li>
<li><a class="reference external" href="#prime-counting-functions">Prime counting functions</a><ul>
<li><a class="reference external" href="#primepi"><tt class="docutils literal"><span class="pre">primepi()</span></tt></a></li>
<li><a class="reference external" href="#primepi2"><tt class="docutils literal"><span class="pre">primepi2()</span></tt></a></li>
<li><a class="reference external" href="#riemannr"><tt class="docutils literal"><span class="pre">riemannr()</span></tt></a></li>
</ul>
</li>
<li><a class="reference external" href="#cyclotomic-polynomials">Cyclotomic polynomials</a><ul>
<li><a class="reference external" href="#cyclotomic"><tt class="docutils literal"><span class="pre">cyclotomic()</span></tt></a></li>
</ul>
</li>
<li><a class="reference external" href="#arithmetic-functions">Arithmetic functions</a><ul>
<li><a class="reference external" href="#mangoldt"><tt class="docutils literal"><span class="pre">mangoldt()</span></tt></a></li>
</ul>
</li>
</ul>
</li>
</ul>

            <h4>Previous topic</h4>
            <p class="topless"><a href="zeta.html"
                                  title="previous chapter">Zeta functions, L-series and polylogarithms</a></p>
            <h4>Next topic</h4>
            <p class="topless"><a href="qfunctions.html"
                                  title="next chapter">q-functions</a></p>
            <h3>This Page</h3>
            <ul class="this-page-menu">
              <li><a href="../_sources/functions/numtheory.txt"
                     rel="nofollow">Show Source</a></li>
            </ul>
          <div id="searchbox" style="display: none">
            <h3>Quick search</h3>
              <form class="search" action="../search.html" method="get">
                <input type="text" name="q" size="18" />
                <input type="submit" value="Go" />
                <input type="hidden" name="check_keywords" value="yes" />
                <input type="hidden" name="area" value="default" />
              </form>
              <p class="searchtip" style="font-size: 90%">
              Enter search terms or a module, class or function name.
              </p>
          </div>
          <script type="text/javascript">$('#searchbox').show(0);</script>
        </div>
      </div>
      <div class="clearer"></div>
    </div>
    <div class="related">
      <h3>Navigation</h3>
      <ul>
        <li class="right" style="margin-right: 10px">
          <a href="../genindex.html" title="General Index"
             >index</a></li>
        <li class="right" >
          <a href="../modindex.html" title="Global Module Index"
             >modules</a> |</li>
        <li class="right" >
          <a href="qfunctions.html" title="q-functions"
             >next</a> |</li>
        <li class="right" >
          <a href="zeta.html" title="Zeta functions, L-series and polylogarithms"
             >previous</a> |</li>
        <li><a href="../index.html">mpmath v0.17 documentation</a> &raquo;</li>
          <li><a href="index.html" >Mathematical functions</a> &raquo;</li> 
      </ul>
    </div>
    <div class="footer">
      &copy; Copyright 2010, Fredrik Johansson.
      Last updated on Feb 06, 2011.
      Created using <a href="http://sphinx.pocoo.org/">Sphinx</a> 0.6.6.
    </div>
  </body>
</html>