Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is to exceed 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - Ser.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.1 lcov report (development 30702-bddb8d6928) Lines: 162 164 98.8 %
Date: 2026-02-23 02:23:56 Functions: 21 21 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000, 2012  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation; either version 2 of the License, or (at your option) any later
       8             : version. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : #include "pari.h"
      16             : #include "paripriv.h"
      17             : /*******************************************************************/
      18             : /*                                                                 */
      19             : /*                     Conversion --> t_SER                        */
      20             : /*                                                                 */
      21             : /*******************************************************************/
      22             : /* x RgX in variable t, return x * (1 + O(t^(l-2))) assuming l >= 2 and
      23             :  * e = v_t(x). Length and valuation changed by normalization, stripping
      24             :  * leading zero terms. Shallow if copy = 0, else gc clean */
      25             : static GEN
      26     1861415 : RgX_to_ser_i(GEN x, long l, long e, int copy)
      27             : {
      28     1861415 :   long i = 2, lx = lg(x), vx = varn(x);
      29             :   GEN y;
      30     1861415 :   if (lx == 2) return zeroser(vx, minss(l - 2, e));
      31     1861331 :   if (l <= 2)
      32             :   {
      33          10 :     if (l == 2 && e != LONG_MAX) return zeroser(vx, e);
      34           0 :     pari_err_BUG("RgX_to_ser (l < 2)");
      35             :   }
      36     1861321 :   y = cgetg(l,t_SER);
      37     1861321 :   if (e == LONG_MAX) { e = 1; lx = 3; } /* e.g. Mod(0,3) * x^0 */
      38     1861305 :   else if (e)
      39             :   {
      40       13643 :     long w = e;
      41      350666 :     while (isrationalzero(gel(x,2))) { x++; w--; }
      42       13643 :     lx -= e;
      43       13643 :     if (w)
      44             :     { /* keep type information, e.g. Mod(0,3) + x */
      45          36 :       GEN z = gel(x,2); /* = 0 */
      46          36 :       if (lx <= 2) gel(y,2) = copy? gcopy(z): z;
      47          31 :       else { x += w; gel(y,2) = gadd(gel(x,2), z); }
      48          36 :       i++;
      49             :     }
      50             :   }
      51     1861321 :   y[1] = evalvarn(vx) | evalvalser(e); /* must come here because of LONG_MAX */
      52     1861321 :   if (lx > l) lx = l;
      53             :   /* 2 <= lx <= l */
      54     1861321 :   if (copy)
      55         167 :     for (; i <lx; i++) gel(y,i) = gcopy(gel(x,i));
      56             :   else
      57     7590886 :     for (; i <lx; i++) gel(y,i) = gel(x,i);
      58     8465107 :   for (     ; i < l; i++) gel(y,i) = gen_0;
      59     1861321 :   return normalizeser(y);
      60             : }
      61             : /* enlarge/truncate t_POL x to a t_SER with lg l */
      62             : GEN
      63     1859303 : RgX_to_ser(GEN x, long l) { return RgX_to_ser_i(x, l, RgX_val(x), 0); }
      64             : GEN
      65        1402 : RgX_to_ser_inexact(GEN x, long l)
      66             : {
      67        1402 :   long i, lx = lg(x);
      68        1402 :   int first = 1;
      69        2434 :   for (i = 2; i < lx && gequal0(gel(x,i)); i++) /* RgX_valrem + normalizeser */
      70        1032 :     if (first && !isexactzero(gel(x,i)))
      71             :     {
      72          10 :       pari_warn(warner,"normalizing a series with 0 leading term");
      73          10 :       first = 0;
      74             :     }
      75        1402 :   return RgX_to_ser_i(x, l, i - 2, 0);
      76             : }
      77             : /* *pd t_POL normalized wrt exact zeros; normalize fully, keeping type
      78             :  * information */
      79             : static long
      80        1358 : RgX_valrem_type(GEN *pd, long *warn)
      81             : {
      82        1358 :   GEN d = *pd, z = gel(d,2);
      83             :   long v;
      84        1358 :   if (!gequal0(z)) return 0;
      85          44 :   *warn = 1;
      86          44 :   if (!signe(d)) { *pd = scalarpol_shallow(z, varn(d)); return degpol(d); }
      87          38 :   v = RgX_valrem_inexact(d, &d);
      88          38 :   if (lg(d) > 2)
      89          38 :     gel(d,2) = gadd(gel(d,2), z);
      90             :   else
      91           0 :     d = scalarpol_shallow(z, varn(d));
      92          38 :   *pd = d; return v;
      93             : }
      94             : static GEN
      95         689 : _rfrac_to_ser(GEN x, long l, long copy)
      96             : {
      97         689 :   GEN a = gel(x,1), d = gel(x,2);
      98         689 :   long warn = 0, v = varn(d), e;
      99         689 :   if (l == 2) return zeroser(v, gvaluation(x, pol_x(v)));
     100         679 :   e = - RgX_valrem(d, &d);
     101         679 :   e -= RgX_valrem_type(&d, &warn);
     102         679 :   if (!signe(d)) pari_err_INV("rfrac_to_ser", gel(x,2));
     103         679 :   if (typ(a) != t_POL || varn(a) != v)
     104             :   {
     105         488 :     a = RgX_Rg_mul(RgXn_inv(d, l - 2), a);
     106         488 :     e += RgX_valrem_type(&a, &warn);
     107             :   }
     108             :   else
     109             :   {
     110         191 :     e += RgX_valrem(a, &a);
     111         191 :     e += RgX_valrem_type(&a, &warn);
     112         191 :     a = RgXn_div(a, d, l - 2);
     113             :   }
     114         679 :   if (warn) pari_warn(warner,"normalizing a series with 0 leading term");
     115         679 :   a = RgX_to_ser_i(a, l, 0, copy);
     116         679 :   setvalser(a, valser(a) + e); return a;
     117             : }
     118             : GEN
     119          26 : rfrac_to_ser(GEN x, long l) { return _rfrac_to_ser(x, l, 1); }
     120             : GEN
     121         663 : rfrac_to_ser_i(GEN x, long l) { return _rfrac_to_ser(x, l, 0); }
     122             : 
     123             : static GEN
     124        5196 : RgV_to_ser_i(GEN x, long v, long l, int copy)
     125             : {
     126        5196 :   long j, lx = minss(lg(x), l-1);
     127             :   GEN y;
     128        5196 :   if (lx == 1) return zeroser(v, l-2);
     129        5191 :   y = cgetg(l, t_SER); y[1] = evalsigne(1)|evalvarn(v)|evalvalser(0);
     130        5191 :   x--;
     131        5191 :   if (copy)
     132      402425 :     for (j = 2; j <= lx; j++) gel(y,j) = gcopy(gel(x,j));
     133             :   else
     134      246222 :     for (j = 2; j <= lx; j++) gel(y,j) = gel(x,j);
     135        5276 :   for (     ; j < l;   j++) gel(y,j) = gen_0;
     136        5191 :   return normalizeser(y);
     137             : }
     138             : GEN
     139        4986 : RgV_to_ser(GEN x, long v, long l) { return RgV_to_ser_i(x, v, l, 0); }
     140             : 
     141             : /* x a t_SER, prec >= 0 */
     142             : GEN
     143        1238 : sertoser(GEN x, long prec)
     144             : {
     145        1238 :   long i, lx = lg(x), l;
     146             :   GEN y;
     147        1238 :   if (lx == 2) return zeroser(varn(x), prec);
     148        1218 :   l = prec+2; lx = minss(lx, l);
     149        1218 :   y = cgetg(l,t_SER); y[1] = x[1];
     150       88878 :   for (i = 2; i < lx; i++) gel(y,i) = gel(x,i);
     151        1677 :   for (     ; i < l;  i++) gel(y,i) = gen_0;
     152        1218 :   return y;
     153             : }
     154             : 
     155             : /* R(1/x) = x^v * n/d, val(n) = val(d) = 0 */
     156             : long
     157         105 : rfracrecip(GEN *pn, GEN *pd)
     158             : {
     159         105 :   long v = degpol(*pd);
     160         105 :   if (typ(*pn) == t_POL && varn(*pn) == varn(*pd))
     161             :   {
     162          50 :     v -= degpol(*pn);
     163          50 :     (void)RgX_valrem(*pn, pn); *pn = RgX_recip(*pn);
     164             :   }
     165         105 :   (void)RgX_valrem(*pd, pd); *pd = RgX_recip(*pd);
     166         105 :   return v;
     167             : }
     168             : 
     169             : /* R(1/x) + O(x^N) */
     170             : GEN
     171          60 : rfracrecip_to_ser_absolute(GEN R, long N)
     172             : {
     173          60 :   GEN n = gel(R,1), d = gel(R,2);
     174          60 :   long v = rfracrecip(&n, &d); /* R(1/x) = x^v * n/d, val(n) = val(d) = 0 */
     175          60 :   if (N <= v) return zeroser(varn(d), N);
     176          60 :   R = rfrac_to_ser_i(mkrfrac(n, d), N-v+2);
     177          60 :   setvalser(R, v); return R;
     178             : }
     179             : 
     180             : /* assume prec >= 0 */
     181             : GEN
     182       26088 : scalarser(GEN x, long v, long prec)
     183             : {
     184             :   long i, l, s;
     185             :   GEN y;
     186             : 
     187       26088 :   if (isexactzero(x))
     188             :   {
     189        1307 :     if (isrationalzero(x)) return zeroser(v, prec);
     190         397 :     y = cgetg(3, t_SER);
     191         397 :     y[1] = evalsigne(0) | _evalvalser(prec) | evalvarn(v);
     192         397 :     gel(y,2) = gcopy(x); return y;
     193             :   }
     194       24781 :   l = prec + 2; y = cgetg(l, t_SER); s = !gequal0(x);
     195       24781 :   y[1] = evalsigne(s) | _evalvalser(0) | evalvarn(v);
     196       61119 :   gel(y,2) = gcopy(x); for (i=3; i<l; i++) gel(y,i) = gen_0;
     197       24781 :   return y;
     198             : }
     199             : 
     200             : /* assume x a t_[SER|POL], apply gtoser to all coeffs */
     201             : static GEN
     202           5 : coefstoser(GEN x, long v, long prec)
     203             : {
     204             :   long i, lx;
     205           5 :   GEN y = cgetg_copy(x, &lx); y[1] = x[1];
     206          15 :   for (i=2; i<lx; i++) gel(y,i) = gtoser(gel(x,i), v, prec);
     207           5 :   return y;
     208             : }
     209             : 
     210             : static void
     211          10 : err_ser_priority(GEN x, long v) { pari_err_PRIORITY("Ser", x, "<", v); }
     212             : /* x a t_POL */
     213             : static GEN
     214          46 : poltoser(GEN x, long v, long prec)
     215             : {
     216          46 :   long s = varncmp(varn(x), v);
     217          46 :   if (s < 0) err_ser_priority(x,v);
     218          41 :   if (s > 0) return scalarser(x, v, prec);
     219          31 :   return RgX_to_ser_i(x, prec+2, RgX_val(x), 1);
     220             : }
     221             : /* x a t_RFRAC */
     222             : static GEN
     223          56 : rfractoser(GEN x, long v, long prec)
     224             : {
     225          56 :   long s = varncmp(varn(gel(x,2)), v);
     226             :   pari_sp av;
     227          56 :   if (s < 0) err_ser_priority(x,v);
     228          51 :   if (s > 0) return scalarser(x, v, prec);
     229          26 :   av = avma; return gc_upto(av, rfrac_to_ser(x, prec+2));
     230             : }
     231             : GEN
     232     3558809 : toser_i(GEN x)
     233             : {
     234     3558809 :   switch(typ(x))
     235             :   {
     236      122378 :     case t_SER: return x;
     237        1402 :     case t_POL: return RgX_to_ser_inexact(x, precdl+2);
     238         118 :     case t_RFRAC: return rfrac_to_ser_i(x, precdl+2);
     239             :   }
     240     3434911 :   return NULL;
     241             : }
     242             : 
     243             : /* conversion: prec ignored if t_VEC or t_SER in variable v */
     244             : GEN
     245         328 : gtoser(GEN x, long v, long prec)
     246             : {
     247         328 :   long tx = typ(x);
     248             : 
     249         328 :   if (v < 0) v = 0;
     250         328 :   if (prec < 0) pari_err_DOMAIN("Ser", "precision", "<", gen_0, stoi(prec));
     251         328 :   if (tx == t_SER)
     252             :   {
     253          25 :     long s = varncmp(varn(x), v);
     254          25 :     if      (s < 0) return coefstoser(x, v, prec);
     255          20 :     else if (s > 0) return scalarser(x, v, prec);
     256          10 :     return gcopy(x);
     257             :   }
     258         303 :   if (is_scalar_t(tx)) return scalarser(x,v,prec);
     259         272 :   switch(tx)
     260             :   {
     261          46 :     case t_POL: return poltoser(x, v, prec);
     262          56 :     case t_RFRAC: return rfractoser(x, v, prec);
     263          30 :     case t_QFB: return RgV_to_ser_i(x, v, 4+1, 1);
     264          10 :     case t_VECSMALL: x = zv_to_ZV(x);/*fall through*/
     265         135 :     case t_VEC: case t_COL:
     266         135 :       if (varncmp(gvar(x), v) <= 0) pari_err_PRIORITY("Ser", x, "<=", v);
     267         130 :       return RgV_to_ser_i(x, v, lg(x)+1, 1);
     268             :   }
     269           5 :   pari_err_TYPE("Ser",x);
     270             :   return NULL; /* LCOV_EXCL_LINE */
     271             : }
     272             : /* impose prec */
     273             : GEN
     274         125 : gtoser_prec(GEN x, long v, long prec)
     275             : {
     276         125 :   pari_sp av = avma;
     277         125 :   if (v < 0) v = 0;
     278         125 :   if (prec < 0) pari_err_DOMAIN("Ser", "precision", "<", gen_0, stoi(prec));
     279         120 :   switch(typ(x))
     280             :   {
     281          20 :     case t_SER: if (varn(x) != v) break;
     282          15 :                 return gc_GEN(av, sertoser(x, prec));
     283          20 :     case t_QFB:
     284          20 :       x = RgV_to_ser_i(mkvec3(gel(x,1),gel(x,2),gel(x,3)), v, prec+2, 1);
     285          20 :       return gc_upto(av, x);
     286          10 :     case t_VECSMALL: x = zv_to_ZV(x);/*fall through*/
     287          30 :     case t_VEC: case t_COL:
     288          30 :       if (varncmp(gvar(x), v) <= 0) pari_err_PRIORITY("Ser", x, "<=", v);
     289          30 :       return RgV_to_ser_i(x, v, prec+2, 1);
     290             :   }
     291          55 :   return gtoser(x, v, prec);
     292             : }
     293             : GEN
     294         376 : Ser0(GEN x, long v, GEN d, long prec)
     295             : {
     296         376 :   if (!d) return gtoser(x, v, prec);
     297         125 :   if (typ(d) != t_INT)
     298             :   {
     299           5 :     d = gceil(d);
     300           5 :     if (typ(d) != t_INT) pari_err_TYPE("Ser [precision]",d);
     301             :   }
     302         125 :   return gtoser_prec(x, v, itos(d));
     303             : }

Generated by: LCOV version 1.16