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-bordeaux1.fr machine (x86_64 architecture), and agregate them in the final report:

The target is 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 - ZG.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 16741-1378b1c) Lines: 31 57 54.4 %
Date: 2014-08-17 Functions: 8 12 66.7 %
Legend: Lines: hit not hit | Branches: + taken - not taken # not executed Branches: 17 44 38.6 %

           Branch data     Line data    Source code
       1                 :            : /* $Id$
       2                 :            : 
       3                 :            : Copyright (C) 2011  The PARI group.
       4                 :            : 
       5                 :            : This file is part of the PARI/GP package.
       6                 :            : 
       7                 :            : PARI/GP is free software; you can redistribute it and/or modify it under the
       8                 :            : terms of the GNU General Public License as published by the Free Software
       9                 :            : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
      10                 :            : ANY WARRANTY WHATSOEVER.
      11                 :            : 
      12                 :            : Check the License for details. You should have received a copy of it, along
      13                 :            : with the package; see the file 'COPYING'. If not, write to the Free Software
      14                 :            : Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. */
      15                 :            : 
      16                 :            : #include "pari.h"
      17                 :            : #include "paripriv.h"
      18                 :            : 
      19                 :            : static int
      20                 :       1100 : cmp_G(void *E, GEN x, GEN y) { (void)E; return cmp_universal(x,y); }
      21                 :            : 
      22                 :            : /* a ZG is either a t_INT or a t_VEC of pairs [g,e] representing
      23                 :            :  * \sum e_i [g_i], e_i in Z, g_i in G. */
      24                 :            : GEN
      25                 :        765 : ZG_normalize(GEN x)
      26                 :            : {
      27         [ -  + ]:        765 :   if (typ(x) == t_INT) return x;
      28                 :        765 :   return sort_factor(x, NULL, &cmp_G);
      29                 :            : }
      30                 :            : GEN
      31                 :        625 : ZG_add(GEN x, GEN y)
      32                 :            : {
      33         [ +  + ]:        625 :   if (typ(x) == t_INT)
      34                 :            :   {
      35         [ +  - ]:        365 :     if (!signe(x)) return y;
      36         [ #  # ]:          0 :     if (typ(y) == t_INT)
      37                 :            :     {
      38         [ #  # ]:          0 :       if (!signe(y)) return x;
      39                 :          0 :       return addii(x,y);
      40                 :            :     }
      41                 :          0 :     x = to_famat_shallow(gen_1,x);
      42                 :            :   }
      43         [ -  + ]:        260 :   else if (typ(y) == t_INT)
      44                 :            :   {
      45         [ #  # ]:          0 :     if (!signe(y)) return x;
      46                 :          0 :     y = to_famat_shallow(gen_1,y);
      47                 :            :   }
      48                 :        260 :   x = merge_factor(x, y, NULL, &cmp_G);
      49         [ +  + ]:        260 :   if (lg(gel(x,1)) == 1) return gen_0;
      50                 :        625 :   return x;
      51                 :            : }
      52                 :            : GEN
      53                 :         75 : ZG_neg(GEN x)
      54                 :            : {
      55         [ -  + ]:         75 :   if (typ(x) == t_INT) return negi(x);
      56                 :         75 :   return mkmat2(gel(x,1),ZC_neg(gel(x,2)));
      57                 :            : }
      58                 :            : GEN
      59                 :         50 : ZG_sub(GEN x, GEN y) { return ZG_add(x, ZG_neg(y)); }
      60                 :            : 
      61                 :            : /* x * c.[1], x in Z[G] */
      62                 :            : GEN
      63                 :          0 : ZG_Z_mul(GEN x, GEN c)
      64                 :            : {
      65 [ #  # ][ #  # ]:          0 :   if (is_pm1(c)) return signe(c) > 0? x: ZG_neg(x);
      66         [ #  # ]:          0 :   if (typ(x) == t_INT) return mulii(x,c);
      67                 :          0 :   return mkmat2(gel(x,1), ZC_Z_mul(gel(x,2), c));
      68                 :            : }
      69                 :            : 
      70                 :            : GEN
      71                 :          0 : ZG_mul(GEN x, GEN y)
      72                 :            : {
      73                 :            :   GEN z, XG, XE;
      74                 :            :   long i, l;
      75         [ #  # ]:          0 :   if (typ(x) == t_INT) return ZG_Z_mul(y, x);
      76         [ #  # ]:          0 :   if (typ(y) == t_INT) return ZG_Z_mul(x, y);
      77                 :          0 :   XG = gel(x,1); XE = gel(x,2); l = lg(XG);
      78                 :          0 :   z = ZG_Z_mul(G_ZG_mul(gel(XG,1), y), gel(XE,1));
      79         [ #  # ]:          0 :   for (i = 2; i < l; i++)
      80                 :          0 :     z = ZG_add(z, ZG_Z_mul(G_ZG_mul(gel(XG,i), y), gel(XE,i)));
      81                 :          0 :   return z;
      82                 :            : }
      83                 :            : #if 0
      84                 :            : static GEN
      85                 :            : ZGV_add(GEN x, GEN y)
      86                 :            : {
      87                 :            :   long i, l;
      88                 :            :   GEN v = cgetg_copy(x, &l);
      89                 :            :   for (i = 1; i < l; i++) gel(v,i) = ZG_add(gel(x,i), gel(y,i));
      90                 :            :   return v;
      91                 :            : }
      92                 :            : static GEN
      93                 :            : ZGV_sub(GEN x, GEN y)
      94                 :            : {
      95                 :            :   long i, l;
      96                 :            :   GEN v = cgetg_copy(x, &l);
      97                 :            :   for (i = 1; i < l; i++) gel(v,i) = ZG_sub(gel(x,i), gel(y,i));
      98                 :            :   return v;
      99                 :            : }
     100                 :            : #endif
     101                 :            : GEN
     102                 :        480 : ZG_G_mul(GEN x, GEN y)
     103                 :            : {
     104                 :            :   long i, l;
     105                 :            :   GEN z, X;
     106         [ +  + ]:        480 :   if (typ(x) == t_INT) return to_famat_shallow(y, x);
     107                 :        290 :   X = gel(x,1);
     108                 :        290 :   z = cgetg_copy(X, &l);
     109         [ +  + ]:        695 :   for (i = 1; i < l; i++) gel(z,i) = gmul(gel(X,i), y);
     110                 :        480 :   return ZG_normalize( mkmat2(z, gel(x,2)) );
     111                 :            : }
     112                 :            : GEN
     113                 :        130 : G_ZG_mul(GEN x, GEN y)
     114                 :            : {
     115                 :            :   long i, l;
     116                 :            :   GEN z, Y;
     117         [ -  + ]:        130 :   if (typ(y) == t_INT) return to_famat_shallow(x, y);
     118                 :        130 :   Y = gel(y,1);
     119                 :        130 :   z = cgetg_copy(Y, &l);
     120         [ +  + ]:        260 :   for (i = 1; i < l; i++) gel(z,i) = gmul(x, gel(Y,i));
     121                 :        130 :   return ZG_normalize( mkmat2(z, gel(y,2)) );
     122                 :            : }
     123                 :            : GEN
     124                 :        210 : ZGC_G_mul(GEN v, GEN x)
     125                 :            : {
     126                 :            :   long i, l;
     127                 :        210 :   GEN w = cgetg_copy(v, &l);
     128         [ +  + ]:        690 :   for (i = 1; i < l; i++) gel(w,i) = ZG_G_mul(gel(v,i), x);
     129                 :        210 :   return w;
     130                 :            : }
     131                 :            : GEN
     132                 :          0 : G_ZGC_mul(GEN x, GEN v)
     133                 :            : {
     134                 :            :   long i, l;
     135                 :          0 :   GEN w = cgetg_copy(v, &l);
     136         [ #  # ]:          0 :   for (i = 1; i < l; i++) gel(w,i) = G_ZG_mul(gel(v,i), x);
     137                 :          0 :   return w;
     138                 :            : }
     139                 :            : GEN
     140                 :          0 : ZGC_Z_mul(GEN v, GEN x)
     141                 :            : {
     142                 :            :   long i, l;
     143                 :          0 :   GEN w = cgetg_copy(v, &l);
     144         [ #  # ]:          0 :   for (i = 1; i < l; i++) gel(w,i) = ZG_Z_mul(gel(v,i), x);
     145                 :          0 :   return w;
     146                 :            : }

Generated by: LCOV version 1.9