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 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 - alglin3.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 20704-afbaa81) Lines: 459 513 89.5 %
Date: 2017-05-26 05:56:59 Functions: 41 48 85.4 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 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. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /********************************************************************/
      15             : /**                                                                **/
      16             : /**                         LINEAR ALGEBRA                         **/
      17             : /**                          (third part)                          **/
      18             : /**                                                                **/
      19             : /********************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : /*******************************************************************/
      24             : /*                                                                 */
      25             : /*                               SUM                               */
      26             : /*                                                                 */
      27             : /*******************************************************************/
      28             : 
      29             : GEN
      30       88637 : vecsum(GEN v)
      31             : {
      32       88637 :   pari_sp av = avma;
      33             :   long i, l;
      34             :   GEN p;
      35       88637 :   if (!is_vec_t(typ(v)))
      36           7 :     pari_err_TYPE("vecsum", v);
      37       88630 :   l = lg(v);
      38       88630 :   if (l == 1) return gen_0;
      39       88623 :   p = gel(v,1);
      40       88623 :   if (l == 2) return gcopy(p);
      41      174909 :   for (i=2; i<l; i++)
      42             :   {
      43      121366 :     p = gadd(p, gel(v,i));
      44      121366 :     if (gc_needed(av, 2))
      45             :     {
      46           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sum");
      47           0 :       p = gerepileupto(av, p);
      48             :     }
      49             :   }
      50       53543 :   return gerepileupto(av, p);
      51             : }
      52             : 
      53             : /*******************************************************************/
      54             : /*                                                                 */
      55             : /*                         TRANSPOSE                               */
      56             : /*                                                                 */
      57             : /*******************************************************************/
      58             : /* A[x0,]~ */
      59             : static GEN
      60     5254008 : row_transpose(GEN A, long x0)
      61             : {
      62     5254008 :   long i, lB = lg(A);
      63     5254008 :   GEN B  = cgetg(lB, t_COL);
      64     5254008 :   for (i=1; i<lB; i++) gel(B, i) = gcoeff(A, x0, i);
      65     5254008 :   return B;
      66             : }
      67             : static GEN
      68       17563 : row_transposecopy(GEN A, long x0)
      69             : {
      70       17563 :   long i, lB = lg(A);
      71       17563 :   GEN B  = cgetg(lB, t_COL);
      72       17563 :   for (i=1; i<lB; i++) gel(B, i) = gcopy(gcoeff(A, x0, i));
      73       17563 :   return B;
      74             : }
      75             : 
      76             : /* No copy*/
      77             : GEN
      78     1023144 : shallowtrans(GEN x)
      79             : {
      80             :   long i, dx, lx;
      81             :   GEN y;
      82     1023144 :   switch(typ(x))
      83             :   {
      84          42 :     case t_VEC: y = leafcopy(x); settyp(y,t_COL); break;
      85        6195 :     case t_COL: y = leafcopy(x); settyp(y,t_VEC); break;
      86             :     case t_MAT:
      87     1016907 :       lx = lg(x); if (lx==1) return cgetg(1,t_MAT);
      88     1016802 :       dx = lgcols(x); y = cgetg(dx,t_MAT);
      89     1016802 :       for (i = 1; i < dx; i++) gel(y,i) = row_transpose(x,i);
      90     1016802 :       break;
      91           0 :     default: pari_err_TYPE("shallowtrans",x); return NULL;
      92             :   }
      93     1023039 :   return y;
      94             : }
      95             : 
      96             : GEN
      97        4949 : gtrans(GEN x)
      98             : {
      99             :   long i, dx, lx;
     100             :   GEN y;
     101        4949 :   switch(typ(x))
     102             :   {
     103         140 :     case t_VEC: y = gcopy(x); settyp(y,t_COL); break;
     104        2905 :     case t_COL: y = gcopy(x); settyp(y,t_VEC); break;
     105             :     case t_MAT:
     106        1897 :       lx = lg(x); if (lx==1) return cgetg(1,t_MAT);
     107        1890 :       dx = lgcols(x); y = cgetg(dx,t_MAT);
     108        1890 :       for (i = 1; i < dx; i++) gel(y,i) = row_transposecopy(x,i);
     109        1890 :       break;
     110           7 :     default: pari_err_TYPE("gtrans",x); return NULL;
     111             :   }
     112        4935 :   return y;
     113             : }
     114             : 
     115             : /*******************************************************************/
     116             : /*                                                                 */
     117             : /*                           EXTRACTION                            */
     118             : /*                                                                 */
     119             : /*******************************************************************/
     120             : 
     121             : static long
     122         182 : str_to_long(char *s, char **pt)
     123             : {
     124         182 :   long a = atol(s);
     125         182 :   while (isspace((int)*s)) s++;
     126         182 :   if (*s == '-' || *s == '+') s++;
     127         182 :   while (isdigit((int)*s) || isspace((int)*s)) s++;
     128         182 :   *pt = s; return a;
     129             : }
     130             : 
     131             : static int
     132         112 : get_range(char *s, long *a, long *b, long *cmpl, long lx)
     133             : {
     134         112 :   long max = lx - 1;
     135             : 
     136         112 :   *a = 1; *b = max;
     137         112 :   if (*s == '^') { *cmpl = 1; s++; } else *cmpl = 0;
     138         112 :   if (!*s) return 0;
     139         112 :   if (*s != '.')
     140             :   {
     141         105 :     *a = str_to_long(s, &s);
     142         105 :     if (*a < 0) *a += lx;
     143         105 :     if (*a<1 || *a>max) return 0;
     144             :   }
     145         112 :   if (*s == '.')
     146             :   {
     147         105 :     s++; if (*s != '.') return 0;
     148         105 :     do s++; while (isspace((int)*s));
     149         105 :     if (*s)
     150             :     {
     151          77 :       *b = str_to_long(s, &s);
     152          77 :       if (*b < 0) *b += lx;
     153          77 :       if (*b<1 || *b>max || *s) return 0;
     154             :     }
     155          98 :     return 1;
     156             :   }
     157           7 :   if (*s) return 0;
     158           7 :   *b = *a; return 1;
     159             : }
     160             : 
     161             : static int
     162          35 : extract_selector_ok(long lx, GEN L)
     163             : {
     164             :   long i, l;
     165          35 :   switch (typ(L))
     166             :   {
     167             :     case t_INT: {
     168             :       long maxj;
     169           7 :       if (!signe(L)) return 1;
     170           7 :       l = lgefint(L)-1;
     171           7 :       maxj = BITS_IN_LONG - bfffo(*int_MSW(L));
     172           7 :       return ((l-2) * BITS_IN_LONG + maxj < lx);
     173             :     }
     174             :     case t_STR: {
     175             :       long first, last, cmpl;
     176           7 :       return get_range(GSTR(L), &first, &last, &cmpl, lx);
     177             :     }
     178             :     case t_VEC: case t_COL:
     179          14 :       l = lg(L);
     180          28 :       for (i=1; i<l; i++)
     181             :       {
     182          21 :         long j = itos(gel(L,i));
     183          21 :         if (j>=lx || j<=0) return 0;
     184             :       }
     185           7 :       return 1;
     186             :     case t_VECSMALL:
     187           7 :       l = lg(L);
     188          21 :       for (i=1; i<l; i++)
     189             :       {
     190          14 :         long j = L[i];
     191          14 :         if (j>=lx || j<=0) return 0;
     192             :       }
     193           7 :       return 1;
     194             :   }
     195           0 :   return 0;
     196             : }
     197             : 
     198             : GEN
     199           0 : shallowmatextract(GEN x, GEN l1, GEN l2)
     200             : {
     201           0 :   long i, j, n1 = lg(l1), n2 = lg(l2);
     202           0 :   GEN M = cgetg(n2, t_MAT);
     203           0 :   for(i=1; i < n2; i++)
     204             :   {
     205           0 :     long ii = l2[i];
     206           0 :     GEN C = cgetg(n1, t_COL);
     207           0 :     for (j=1; j < n1; j++)
     208             :     {
     209           0 :       long jj = l1[j];
     210           0 :       gel(C, j) = gmael(x, ii, jj);
     211             :     }
     212           0 :     gel(M, i) = C;
     213             :   }
     214           0 :   return M;
     215             : }
     216             : 
     217             : GEN
     218       21602 : shallowextract(GEN x, GEN L)
     219             : {
     220       21602 :   long i,j, tl = typ(L), tx = typ(x), lx = lg(x);
     221             :   GEN y;
     222             : 
     223       21602 :   switch(tx)
     224             :   {
     225             :     case t_VEC:
     226             :     case t_COL:
     227             :     case t_MAT:
     228       21595 :     case t_VECSMALL: break;
     229           7 :     default: pari_err_TYPE("extract",x);
     230             : 
     231             :   }
     232       21595 :   if (tl==t_INT)
     233             :   { /* extract components of x as per the bits of mask L */
     234             :     long k, l, ix, iy, maxj;
     235             :     GEN Ld;
     236        3066 :     if (!signe(L)) return cgetg(1,tx);
     237        3059 :     y = new_chunk(lx);
     238        3059 :     l = lgefint(L)-1; ix = iy = 1;
     239        3059 :     maxj = BITS_IN_LONG - bfffo(*int_MSW(L));
     240        3059 :     if ((l-2) * BITS_IN_LONG + maxj >= lx)
     241           7 :       pari_err_TYPE("vecextract [mask too large]", L);
     242        3427 :     for (k = 2, Ld = int_LSW(L); k < l; k++, Ld = int_nextW(Ld))
     243             :     {
     244         375 :       ulong B = *Ld;
     245       20439 :       for (j = 0; j < BITS_IN_LONG; j++, B >>= 1, ix++)
     246       20064 :         if (B & 1) y[iy++] = x[ix];
     247             :     }
     248             :     { /* k = l */
     249        3052 :       ulong B = *Ld;
     250       27319 :       for (j = 0; j < maxj; j++, B >>= 1, ix++)
     251       24267 :         if (B & 1) y[iy++] = x[ix];
     252             :     }
     253        3052 :     y[0] = evaltyp(tx) | evallg(iy);
     254        3052 :     return y;
     255             :   }
     256       18529 :   if (tl==t_STR)
     257             :   {
     258         105 :     char *s = GSTR(L);
     259             :     long first, last, cmpl, d;
     260         105 :     if (! get_range(s, &first, &last, &cmpl, lx))
     261           7 :       pari_err_TYPE("vecextract [incorrect range]", L);
     262          98 :     if (lx == 1) return cgetg(1,tx);
     263          98 :     d = last - first;
     264          98 :     if (cmpl)
     265             :     {
     266          21 :       if (d >= 0)
     267             :       {
     268          14 :         y = cgetg(lx - (1+d),tx);
     269          14 :         for (j=1; j<first; j++) gel(y,j) = gel(x,j);
     270          14 :         for (i=last+1; i<lx; i++,j++) gel(y,j) = gel(x,i);
     271             :       }
     272             :       else
     273             :       {
     274           7 :         y = cgetg(lx - (1-d),tx);
     275           7 :         for (j=1,i=lx-1; i>first; i--,j++) gel(y,j) = gel(x,i);
     276           7 :         for (i=last-1; i>0; i--,j++) gel(y,j) = gel(x,i);
     277             :       }
     278             :     }
     279             :     else
     280             :     {
     281          77 :       if (d >= 0)
     282             :       {
     283          35 :         y = cgetg(d+2,tx);
     284          35 :         for (i=first,j=1; i<=last; i++,j++) gel(y,j) = gel(x,i);
     285             :       }
     286             :       else
     287             :       {
     288          42 :         y = cgetg(2-d,tx);
     289          42 :         for (i=first,j=1; i>=last; i--,j++) gel(y,j) = gel(x,i);
     290             :       }
     291             :     }
     292          98 :     return y;
     293             :   }
     294             : 
     295       18424 :   if (is_vec_t(tl))
     296             :   {
     297          63 :     long ll=lg(L); y=cgetg(ll,tx);
     298         140 :     for (i=1; i<ll; i++)
     299             :     {
     300          91 :       j = itos(gel(L,i));
     301          91 :       if (j<=0) pari_err_COMPONENT("vecextract","<=",gen_0,stoi(j));
     302          84 :       if (j>=lx) pari_err_COMPONENT("vecextract",">=",stoi(lx),stoi(j));
     303          77 :       gel(y,i) = gel(x,j);
     304             :     }
     305          49 :     return y;
     306             :   }
     307       18361 :   if (tl == t_VECSMALL)
     308             :   {
     309       18354 :     long ll=lg(L); y=cgetg(ll,tx);
     310       76818 :     for (i=1; i<ll; i++)
     311             :     {
     312       58464 :       j = L[i];
     313       58464 :       if (j<=0) pari_err_COMPONENT("vecextract","<=",gen_0,stoi(j));
     314       58464 :       if (j>=lx) pari_err_COMPONENT("vecextract",">=",stoi(lx),stoi(j));
     315       58464 :       gel(y,i) = gel(x,j);
     316             :     }
     317       18354 :     return y;
     318             :   }
     319           7 :   pari_err_TYPE("vecextract [mask]", L);
     320             :   return NULL; /* LCOV_EXCL_LINE */
     321             : }
     322             : 
     323             : /* does the component selector l select 0 component ? */
     324             : static int
     325          77 : select_0(GEN l)
     326             : {
     327          77 :   switch(typ(l))
     328             :   {
     329             :     case t_INT:
     330          14 :       return (!signe(l));
     331             :     case t_VEC: case t_COL: case t_VECSMALL:
     332          42 :       return (lg(l) == 1);
     333             :   }
     334          21 :   return 0;
     335             : }
     336             : 
     337             : GEN
     338       18921 : extract0(GEN x, GEN l1, GEN l2)
     339             : {
     340       18921 :   pari_sp av = avma, av2;
     341             :   GEN y;
     342       18921 :   if (! l2)
     343             :   {
     344       18844 :     y = shallowextract(x, l1);
     345       18802 :     if (lg(y) == 1 || typ(y) == t_VECSMALL) return y;
     346       18795 :     av2 = avma;
     347       18795 :     y = gcopy(y);
     348             :   }
     349             :   else
     350             :   {
     351          77 :     if (typ(x) != t_MAT) pari_err_TYPE("extract",x);
     352          77 :     y = shallowextract(x,l2);
     353          77 :     if (select_0(l1)) { avma = av; return zeromat(0, lg(y)-1); }
     354          63 :     if (lg(y) == 1 && lg(x) > 1)
     355             :     {
     356          35 :       if (!extract_selector_ok(lgcols(x), l1))
     357           7 :         pari_err_TYPE("vecextract [incorrect mask]", l1);
     358          28 :       avma = av; return cgetg(1, t_MAT);
     359             :     }
     360          28 :     y = shallowextract(shallowtrans(y), l1);
     361          28 :     av2 = avma;
     362          28 :     y = gtrans(y);
     363             :   }
     364       18823 :   stackdummy(av, av2);
     365       18823 :   return y;
     366             : }
     367             : 
     368             : static long
     369        3662 : vecslice_parse_arg(long lA, long *y1, long *y2, long *skip)
     370             : {
     371        3662 :   *skip=0;
     372        3662 :   if (*y1==LONG_MAX)
     373             :   {
     374         119 :     if (*y2!=LONG_MAX)
     375             :     {
     376          63 :       if (*y2<0) *y2 += lA;
     377          63 :       if (*y2<0 || *y2==LONG_MAX || *y2>=lA)
     378           0 :         pari_err_DIM("_[..]");
     379          63 :       *skip=*y2;
     380             :     }
     381         119 :     *y1 = 1; *y2 = lA-1;
     382             :   }
     383        3543 :   else if (*y2==LONG_MAX) *y2 = *y1;
     384        3662 :   if (*y1<=0) *y1 += lA;
     385        3662 :   if (*y2<0) *y2 += lA;
     386        3662 :   if (*y1<=0 || *y1>*y2+1 || *y2>=lA) pari_err_DIM("_[..]");
     387        3648 :   return *y2 - *y1 + 2 - !!*skip;
     388             : }
     389             : 
     390             : static GEN
     391        3844 : vecslice_i(GEN A, long t, long lB, long y1, long skip)
     392             : {
     393        3844 :   GEN B = cgetg(lB, t);
     394             :   long i;
     395      108464 :   for (i=1; i<lB; i++, y1++)
     396             :   {
     397      104620 :     if (y1 == skip) { i--; continue; }
     398      104557 :     gel(B,i) = gcopy(gel(A,y1));
     399             :   }
     400        3844 :   return B;
     401             : }
     402             : 
     403             : static GEN
     404          14 : rowslice_i(GEN A, long lB, long x1, long y1, long skip)
     405             : {
     406          14 :   GEN B = cgetg(lB, t_VEC);
     407             :   long i;
     408          77 :   for (i=1; i<lB; i++, y1++)
     409             :   {
     410          63 :     if (y1 == skip) { i--; continue; }
     411          56 :     gel(B,i) = gcopy(gcoeff(A,x1,y1));
     412             :   }
     413          14 :   return B;
     414             : }
     415             : 
     416             : static GEN
     417           0 : rowsmallslice_i(GEN A, long lB, long x1, long y1, long skip)
     418             : {
     419           0 :   GEN B = cgetg(lB, t_VECSMALL);
     420             :   long i;
     421           0 :   for (i=1; i<lB; i++, y1++)
     422             :   {
     423           0 :     if (y1 == skip) { i--; continue; }
     424           0 :     B[i] = coeff(A,x1,y1);
     425             :   }
     426           0 :   return B;
     427             : }
     428             : 
     429             : static GEN
     430          28 : vecsmallslice_i(GEN A, long t, long lB, long y1, long skip)
     431             : {
     432          28 :   GEN B = cgetg(lB, t);
     433             :   long i;
     434         126 :   for (i=1; i<lB; i++, y1++)
     435             :   {
     436          98 :     if (y1 == skip) { i--; continue; }
     437          91 :     B[i] = A[y1];
     438             :   }
     439          28 :   return B;
     440             : }
     441             : GEN
     442        3382 : vecslice0(GEN A, long y1, long y2)
     443             : {
     444        3382 :   long skip, lB, t = typ(A);
     445        3382 :   lB = vecslice_parse_arg(lg(A), &y1, &y2, &skip);
     446        3368 :   switch(t)
     447             :   {
     448             :     case t_VEC: case t_COL:
     449        3333 :       return vecslice_i(A, t,lB,y1,skip);
     450             :     case t_VECSMALL:
     451          28 :       return vecsmallslice_i(A, t,lB,y1,skip);
     452             :     default:
     453           7 :       pari_err_TYPE("_[_.._]",A);
     454           0 :       return NULL;
     455             :   }
     456             : }
     457             : 
     458             : GEN
     459         147 : matslice0(GEN A, long x1, long x2, long y1, long y2)
     460             : {
     461             :   GEN B;
     462         147 :   long i, lB, lA = lg(A), rA, t, skip, rskip, rlB;
     463         147 :   long is_col = y1!=LONG_MAX && y2==LONG_MAX;
     464         147 :   long is_row = x1!=LONG_MAX && x2==LONG_MAX;
     465             :   GEN (*slice)(GEN A, long t, long lB, long y1, long skip);
     466         147 :   if (typ(A)!=t_MAT) pari_err_TYPE("_[_.._,_.._]",A);
     467         147 :   lB = vecslice_parse_arg(lA, &y1, &y2, &skip);
     468         147 :   if (is_col) return vecslice0(gel(A, y1), x1, x2);
     469         133 :   rA = lg(A)==1 ? 1: lgcols(A);
     470         133 :   rlB = vecslice_parse_arg(rA, &x1, &x2, &rskip);
     471         133 :   t = lg(A)==1 ? t_COL: typ(gel(A,1));
     472         133 :   if (is_row) return t == t_COL ? rowslice_i(A, lB, x1, y1, skip):
     473           0 :                                   rowsmallslice_i(A, lB, x1, y1, skip);
     474         119 :   slice = t == t_COL? &vecslice_i: &vecsmallslice_i;
     475             : 
     476         119 :   B = cgetg(lB, t_MAT);
     477         644 :   for (i=1; i<lB; i++, y1++)
     478             :   {
     479         525 :     if (y1 == skip) { i--; continue; }
     480         511 :     gel(B,i) = slice(gel(A,y1),t,rlB, x1, rskip);
     481             :   }
     482         119 :   return B;
     483             : }
     484             : 
     485             : GEN
     486        1129 : vecrange(GEN a, GEN b)
     487             : {
     488             :   GEN y;
     489             :   long i, l;
     490        1129 :   if (typ(a)!=t_INT) pari_err_TYPE("[_.._]",a);
     491        1122 :   if (typ(b)!=t_INT) pari_err_TYPE("[_.._]",b);
     492        1115 :   if (cmpii(a,b)>0) return cgetg(1,t_VEC);
     493        1108 :   l = itos(subii(b,a))+1;
     494        1108 :   a = setloop(a);
     495        1108 :   y = cgetg(l+1, t_VEC);
     496      260213 :   for (i=1; i<=l; a = incloop(a), i++)
     497      259105 :     gel(y,i) = icopy(a);
     498        1108 :   return y;
     499             : }
     500             : 
     501             : GEN
     502           0 : vecrangess(long a, long b)
     503             : {
     504             :   GEN y;
     505             :   long i, l;
     506           0 :   if (a>b) return cgetg(1,t_VEC);
     507           0 :   l = b-a+1;
     508           0 :   y = cgetg(l+1, t_VEC);
     509           0 :   for (i=1; i<=l; a++, i++)
     510           0 :     gel(y,i) = stoi(a);
     511           0 :   return y;
     512             : }
     513             : 
     514             : GEN
     515          84 : genindexselect(void *E, long (*f)(void* E, GEN x), GEN A)
     516             : {
     517             :   long l, i, lv;
     518             :   GEN v, z;
     519             :   pari_sp av;
     520          84 :   clone_lock(A);
     521          84 :   switch(typ(A))
     522             :   {
     523             :     case t_LIST:
     524           7 :       z = list_data(A);
     525           7 :       l = z? lg(z): 1;
     526           7 :       break;
     527             :     case t_VEC: case t_COL: case t_MAT:
     528          70 :       l = lg(A);
     529          70 :       z = A;
     530          70 :       break;
     531             :     default:
     532           7 :       pari_err_TYPE("select",A);
     533             :       return NULL;/*LCOV_EXCL_LINE*/
     534             :   }
     535          77 :   v = cgetg(l, t_VECSMALL);
     536          77 :   av = avma;
     537         749 :   for (i = lv = 1; i < l; i++) {
     538         672 :     if (f(E, gel(z,i))) v[lv++] = i;
     539         672 :     avma = av;
     540             :   }
     541          77 :   clone_unlock(A); fixlg(v, lv); return v;
     542             : }
     543             : static GEN
     544          76 : extract_copy(GEN A, GEN v)
     545             : {
     546          76 :   long i, l = lg(v);
     547          76 :   GEN B = cgetg(l, typ(A));
     548          76 :   for (i = 1; i < l; i++) gel(B,i) = gcopy(gel(A,v[i]));
     549          76 :   return B;
     550             : }
     551             : /* as genselect, but treat A [ t_VEC,t_COL, or t_MAT] as a t_VEC */
     552             : GEN
     553           0 : vecselect(void *E, long (*f)(void* E, GEN x), GEN A)
     554             : {
     555             :   GEN v;
     556           0 :   clone_lock(A);
     557           0 :   v = genindexselect(E, f, A);
     558           0 :   A = extract_copy(A, v); settyp(A, t_VEC);
     559           0 :   clone_unlock(A); return A;
     560             : }
     561             : GEN
     562          77 : genselect(void *E, long (*f)(void* E, GEN x), GEN A)
     563             : {
     564             :   GEN y, z, v;/* v left on stack for efficiency */
     565          77 :   clone_lock(A);
     566          77 :   switch(typ(A))
     567             :   {
     568             :     case t_LIST:
     569          14 :       z = list_data(A);
     570          14 :       if (!z) y = mklist();
     571             :       else
     572             :       {
     573             :         GEN B;
     574          14 :         y = cgetg(3, t_LIST);
     575          14 :         v = genindexselect(E, f, z);
     576          14 :         B = extract_copy(z, v);
     577          14 :         y[1] = lg(B)-1;
     578          14 :         list_data(y) = B;
     579             :       }
     580          14 :       break;
     581             :     case t_VEC: case t_COL: case t_MAT:
     582          56 :       v = genindexselect(E, f, A);
     583          56 :       y = extract_copy(A, v);
     584          56 :       break;
     585             :     default:
     586           7 :       pari_err_TYPE("select",A);
     587             :       return NULL;/*LCOV_EXCL_LINE*/
     588             :   }
     589          70 :   clone_unlock(A); return y;
     590             : }
     591             : 
     592             : static void
     593        8488 : check_callgen1(GEN f, const char *s)
     594             : {
     595        8488 :   if (typ(f) != t_CLOSURE || closure_is_variadic(f)  || closure_arity(f) < 1)
     596           0 :     pari_err_TYPE(s, f);
     597        8488 : }
     598             : 
     599             : GEN
     600          91 : select0(GEN f, GEN x, long flag)
     601             : {
     602          91 :   check_callgen1(f, "select");
     603          91 :   switch(flag)
     604             :   {
     605          77 :     case 0: return genselect((void *) f, gp_callbool, x);
     606          14 :     case 1: return genindexselect((void *) f, gp_callbool, x);
     607           0 :     default: pari_err_FLAG("select");
     608             :              return NULL;/*LCOV_EXCL_LINE*/
     609             :   }
     610             : }
     611             : 
     612             : GEN
     613          12 : parselect(GEN C, GEN D, long flag)
     614             : {
     615             :   pari_sp av, av2;
     616          12 :   long lv, l = lg(D), i, pending = 0, workid;
     617             :   GEN V, worker, done;
     618             :   struct pari_mt pt;
     619          12 :   check_callgen1(C, "parselect");
     620          12 :   if (!is_vec_t(typ(D))) pari_err_TYPE("parapply",D);
     621          12 :   V = cgetg(l, t_VECSMALL); av = avma;
     622          12 :   worker = strtoclosure("_parapply_worker", 1, C);
     623          12 :   av2 = avma;
     624          12 :   mt_queue_start_lim(&pt, worker, l-1);
     625       12084 :   for (i=1; i<l || pending; i++)
     626             :   {
     627       12072 :     mt_queue_submit(&pt, i, i<l? mkvec(gel(D,i)): NULL);
     628       12072 :     done = mt_queue_get(&pt, &workid, &pending);
     629       12072 :     if (done) V[workid] = !gequal0(done);
     630       12072 :     avma = av2;
     631             :   }
     632          12 :   mt_queue_end(&pt);
     633          12 :   avma = av;
     634       12024 :   for (lv=1, i=1; i<l; i++)
     635       12012 :     if (V[i]) V[lv++]=i;
     636          12 :   fixlg(V, lv);
     637          12 :   return flag? V: extract_copy(D, V);
     638             : }
     639             : 
     640             : GEN
     641           0 : veccatapply(void *E, GEN (*f)(void* E, GEN x), GEN x)
     642             : {
     643           0 :   pari_sp av = avma;
     644           0 :   GEN v = vecapply(E, f, x);
     645           0 :   return lg(v) == 1? v: gerepilecopy(av, shallowconcat1(v));
     646             : }
     647             : 
     648             : static GEN
     649          14 : vecapply2(void *E, GEN (*f)(void* E, GEN x), GEN x)
     650             : {
     651             :   long i, lx;
     652          14 :   GEN y = cgetg_copy(x, &lx); y[1] = x[1];
     653          14 :   for (i=2; i<lx; i++) gel(y,i) = f(E, gel(x,i));
     654          14 :   return y;
     655             : }
     656             : static GEN
     657      114163 : vecapply1(void *E, GEN (*f)(void* E, GEN x), GEN x)
     658             : {
     659             :   long i, lx;
     660      114163 :   GEN y = cgetg_copy(x, &lx);
     661      114163 :   for (i=1; i<lx; i++) gel(y,i) = f(E, gel(x,i));
     662      114163 :   return y;
     663             : }
     664             : static GEN
     665           7 : mapapply1(void *E, GEN (*f)(void* E, GEN x), GEN x)
     666             : {
     667             :   long i, lx;
     668           7 :   GEN y = cgetg_copy(x, &lx);
     669          28 :   for (i=1; i<lx; i++)
     670             :   {
     671          21 :     GEN xi = gel(x, i);
     672          42 :     gel(y,i) = mkvec2(mkvec2(gcopy(gmael(xi, 1, 1)), f(E,gmael(xi, 1, 2))),
     673          21 :                              gcopy(gel(xi, 2)));
     674             :   }
     675           7 :   return y;
     676             : }
     677             : /* as genapply, but treat A [ t_VEC,t_COL, or t_MAT] as a t_VEC */
     678             : GEN
     679      105770 : vecapply(void *E, GEN (*f)(void* E, GEN x), GEN x)
     680             : {
     681             :   GEN y;
     682      105770 :   clone_lock(x); y = vecapply1(E,f,x);
     683      105770 :   clone_unlock(x); settyp(y, t_VEC); return y;
     684             : }
     685             : GEN
     686        8379 : genapply(void *E, GEN (*f)(void* E, GEN x), GEN x)
     687             : {
     688        8379 :   long i, lx, tx = typ(x);
     689             :   GEN y, z;
     690        8379 :   if (is_scalar_t(tx)) return f(E, x);
     691        8379 :   clone_lock(x);
     692        8379 :   switch(tx) {
     693           7 :     case t_POL: y = normalizepol(vecapply2(E,f,x)); break;
     694             :     case t_SER:
     695           7 :       y = ser_isexactzero(x)? gcopy(x): normalize(vecapply2(E,f,x));
     696           7 :       break;
     697             :     case t_LIST:
     698             :       {
     699          28 :         long t = list_typ(x);
     700          28 :         z = list_data(x);
     701          28 :         if (!z)
     702           7 :           y = mklist_typ(t);
     703             :         else
     704             :         {
     705          21 :           y = cgetg(3, t_LIST);
     706          21 :           y[1] = evaltyp(t)|evallg(lg(z)-1);
     707          21 :           switch(t)
     708             :           {
     709             :           case t_LIST_RAW:
     710          14 :             list_data(y) = vecapply1(E,f,z);
     711          14 :             break;
     712             :           case t_LIST_MAP:
     713           7 :             list_data(y) = mapapply1(E,f,z);
     714           7 :             break;
     715             :           }
     716             :         }
     717             :       }
     718          28 :       break;
     719             :     case t_MAT:
     720          42 :       y = cgetg_copy(x, &lx);
     721          42 :       for (i = 1; i < lx; i++) gel(y,i) = vecapply1(E,f,gel(x,i));
     722          42 :       break;
     723             : 
     724        8295 :     case t_VEC: case t_COL: y = vecapply1(E,f,x); break;
     725             :     default:
     726           0 :       pari_err_TYPE("apply",x);
     727             :       return NULL;/*LCOV_EXCL_LINE*/
     728             :   }
     729        8379 :   clone_unlock(x); return y;
     730             : }
     731             : 
     732             : GEN
     733        8379 : apply0(GEN f, GEN x)
     734             : {
     735        8379 :   check_callgen1(f, "apply");
     736        8379 :   return genapply((void *) f, gp_call, x);
     737             : }
     738             : 
     739             : GEN
     740         252 : vecselapply(void *Epred, long (*pred)(void* E, GEN x), void *Efun,
     741             :                          GEN (*fun)(void* E, GEN x), GEN A)
     742             : {
     743             :   GEN y;
     744         252 :   long i, l = lg(A), nb=1;
     745         252 :   clone_lock(A); y = cgetg(l, t_VEC);
     746      167139 :   for (i=1; i<l; i++)
     747      166887 :     if (pred(Epred, gel(A,i))) gel(y,nb++) = fun(Efun, gel(A,i));
     748         252 :   fixlg(y,nb); clone_unlock(A); return y;
     749             : }
     750             : 
     751             : GEN
     752           0 : veccatselapply(void *Epred, long (*pred)(void* E, GEN x), void *Efun,
     753             :                             GEN (*fun)(void* E, GEN x), GEN A)
     754             : {
     755           0 :   pari_sp av = avma;
     756           0 :   GEN v = vecselapply(Epred, pred, Efun, fun, A);
     757           0 :   return lg(v) == 1? v: gerepilecopy(av, shallowconcat1(v));
     758             : }
     759             : 
     760             : GEN
     761       12003 : parapply_worker(GEN d, GEN C)
     762             : {
     763       12003 :   return closure_callgen1(C, d);
     764             : }
     765             : 
     766             : GEN
     767           6 : parapply(GEN C, GEN D)
     768             : {
     769           6 :   pari_sp av = avma;
     770           6 :   long l = lg(D), i, pending = 0, workid;
     771             :   GEN V, worker, done;
     772             :   struct pari_mt pt;
     773           6 :   check_callgen1(C, "parapply");
     774           6 :   if (!is_vec_t(typ(D))) pari_err_TYPE("parapply",D);
     775           6 :   worker = strtoclosure("_parapply_worker", 1, C);
     776           6 :   V = cgetg(l, typ(D));
     777           6 :   mt_queue_start_lim(&pt, worker, l-1);
     778          20 :   for (i=1; i<l || pending; i++)
     779             :   {
     780          14 :     mt_queue_submit(&pt, i, i<l? mkvec(gel(D,i)): NULL);
     781          14 :     done = mt_queue_get(&pt, &workid, &pending);
     782          14 :     if (done) gel(V,workid) = done;
     783             :   }
     784           6 :   mt_queue_end(&pt);
     785           6 :   return gerepilecopy(av, V);
     786             : }
     787             : 
     788             : GEN
     789          28 : genfold(void *E, GEN (*f)(void* E, GEN x, GEN y), GEN x)
     790             : {
     791          28 :   pari_sp av = avma;
     792             :   GEN z;
     793          28 :   long i, l = lg(x);
     794          28 :   if (!is_vec_t(typ(x))|| l==1  ) pari_err_TYPE("fold",x);
     795          28 :   clone_lock(x);
     796          28 :   z = gel(x,1);
     797         119 :   for (i=2; i<l; i++)
     798             :   {
     799          91 :     z = f(E,z,gel(x,i));
     800          91 :     if (gc_needed(av, 2))
     801             :     {
     802           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"fold");
     803           0 :       z = gerepilecopy(av, z);
     804             :     }
     805             :   }
     806          28 :   clone_unlock(x);
     807          28 :   return gerepilecopy(av, z);
     808             : }
     809             : 
     810             : GEN
     811          28 : fold0(GEN f, GEN x)
     812             : {
     813          28 :   if (typ(f) != t_CLOSURE || closure_arity(f) < 2) pari_err_TYPE("apply",f);
     814          28 :   return genfold((void *) f, gp_call2, x);
     815             : }
     816             : /*******************************************************************/
     817             : /*                                                                 */
     818             : /*                     SCALAR-MATRIX OPERATIONS                    */
     819             : /*                                                                 */
     820             : /*******************************************************************/
     821             : GEN
     822       53191 : gtomat(GEN x)
     823             : {
     824             :   long lx, i;
     825             :   GEN y;
     826             : 
     827       53191 :   if (!x) return cgetg(1, t_MAT);
     828       53184 :   switch(typ(x))
     829             :   {
     830             :     case t_LIST:
     831          21 :       if (list_typ(x)==t_LIST_MAP)
     832           7 :         return maptomat(x);
     833          14 :       x = list_data(x);
     834          14 :       if (!x) return cgetg(1, t_MAT);
     835             :       /* fall through */
     836             :     case t_VEC: {
     837         294 :       lx=lg(x); y=cgetg(lx,t_MAT);
     838         294 :       if (lx == 1) break;
     839         294 :       if (typ(gel(x,1)) == t_COL) {
     840         133 :         long h = lgcols(x);
     841         868 :         for (i=2; i<lx; i++) {
     842         735 :           if (typ(gel(x,i)) != t_COL || lg(gel(x,i)) != h) break;
     843             :         }
     844         133 :         if (i == lx) { /* matrix with h-1 rows */
     845         133 :           y = cgetg(lx, t_MAT);
     846         133 :           for (i=1 ; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
     847         133 :           return y;
     848             :         }
     849             :       }
     850         161 :       for (i=1; i<lx; i++) gel(y,i) = mkcolcopy(gel(x,i));
     851         161 :       break;
     852             :     }
     853             :     case t_COL:
     854        7194 :       lx = lg(x);
     855        7194 :       if (lx == 1) return cgetg(1, t_MAT);
     856        7194 :       if (typ(gel(x,1)) == t_VEC) {
     857           7 :         long j, h = lg(gel(x,1));
     858          14 :         for (i=2; i<lx; i++) {
     859           7 :           if (typ(gel(x,i)) != t_VEC || lg(gel(x,i)) != h) break;
     860             :         }
     861           7 :         if (i == lx) { /* matrix with h cols */
     862           7 :           y = cgetg(h, t_MAT);
     863          28 :           for (j=1 ; j<h; j++) {
     864          21 :             gel(y,j) = cgetg(lx, t_COL);
     865          21 :             for (i=1; i<lx; i++) gcoeff(y,i,j) = gcopy(gmael(x,i,j));
     866             :           }
     867           7 :           return y;
     868             :         }
     869             :       }
     870        7187 :       y = mkmatcopy(x); break;
     871             :     case t_MAT:
     872       35882 :       y = gcopy(x); break;
     873             :     case t_QFI: case t_QFR: {
     874             :       GEN b;
     875        8428 :       y = cgetg(3,t_MAT); b = gmul2n(gel(x,2),-1);
     876        8428 :       gel(y,1) = mkcol2(icopy(gel(x,1)), b);
     877        8428 :       gel(y,2) = mkcol2(b, icopy(gel(x,3)));
     878        8428 :       break;
     879             :     }
     880             :     default:
     881        1372 :       y = cgetg(2,t_MAT); gel(y,1) = mkcolcopy(x);
     882        1372 :       break;
     883             :   }
     884       53030 :   return y;
     885             : }
     886             : 
     887             : /* create the diagonal matrix, whose diagonal is given by x */
     888             : GEN
     889         525 : diagonal(GEN x)
     890             : {
     891         525 :   long j, lx, tx = typ(x);
     892             :   GEN y;
     893             : 
     894         525 :   if (! is_matvec_t(tx)) return scalarmat(x,1);
     895         518 :   if (tx==t_MAT)
     896             :   {
     897          14 :     if (RgM_isdiagonal(x)) return gcopy(x);
     898           7 :     pari_err_TYPE("diagonal",x);
     899             :   }
     900         504 :   lx=lg(x); y=cgetg(lx,t_MAT);
     901        1295 :   for (j=1; j<lx; j++)
     902             :   {
     903         791 :     gel(y,j) = zerocol(lx-1);
     904         791 :     gcoeff(y,j,j) = gcopy(gel(x,j));
     905             :   }
     906         504 :   return y;
     907             : }
     908             : /* same, assuming x is a t_VEC/t_COL. Not memory clean. */
     909             : GEN
     910       15288 : diagonal_shallow(GEN x)
     911             : {
     912       15288 :   long j, lx = lg(x);
     913       15288 :   GEN y = cgetg(lx,t_MAT);
     914             : 
     915       42322 :   for (j=1; j<lx; j++)
     916             :   {
     917       27034 :     gel(y,j) = zerocol(lx-1);
     918       27034 :     gcoeff(y,j,j) = gel(x,j);
     919             :   }
     920       15288 :   return y;
     921             : }
     922             : 
     923             : /* compute m*diagonal(d) */
     924             : GEN
     925           7 : matmuldiagonal(GEN m, GEN d)
     926             : {
     927             :   long j, lx;
     928           7 :   GEN y = cgetg_copy(m, &lx);
     929             : 
     930           7 :   if (typ(m)!=t_MAT) pari_err_TYPE("matmuldiagonal",m);
     931           7 :   if (! is_vec_t(typ(d))) pari_err_TYPE("matmuldiagonal",d);
     932           7 :   if (lg(d) != lx) pari_err_OP("operation 'matmuldiagonal'", m,d);
     933           7 :   for (j=1; j<lx; j++) gel(y,j) = RgC_Rg_mul(gel(m,j), gel(d,j));
     934           7 :   return y;
     935             : }
     936             : 
     937             : /* compute A*B assuming the result is a diagonal matrix */
     938             : GEN
     939           7 : matmultodiagonal(GEN A, GEN B)
     940             : {
     941           7 :   long i, j, hA, hB, lA = lg(A), lB = lg(B);
     942           7 :   GEN y = matid(lB-1);
     943             : 
     944           7 :   if (typ(A) != t_MAT) pari_err_TYPE("matmultodiagonal",A);
     945           7 :   if (typ(B) != t_MAT) pari_err_TYPE("matmultodiagonal",B);
     946           7 :   hA = (lA == 1)? lB: lgcols(A);
     947           7 :   hB = (lB == 1)? lA: lgcols(B);
     948           7 :   if (lA != hB || lB != hA) pari_err_OP("operation 'matmultodiagonal'", A,B);
     949          56 :   for (i=1; i<lB; i++)
     950             :   {
     951          49 :     GEN z = gen_0;
     952          49 :     for (j=1; j<lA; j++) z = gadd(z, gmul(gcoeff(A,i,j),gcoeff(B,j,i)));
     953          49 :     gcoeff(y,i,i) = z;
     954             :   }
     955           7 :   return y;
     956             : }
     957             : 
     958             : /* [m[1,1], ..., m[l,l]], internal */
     959             : GEN
     960      151536 : RgM_diagonal_shallow(GEN m)
     961             : {
     962      151536 :   long i, lx = lg(m);
     963      151536 :   GEN y = cgetg(lx,t_VEC);
     964      151536 :   for (i=1; i<lx; i++) gel(y, i) = gcoeff(m,i,i);
     965      151536 :   return y;
     966             : }
     967             : 
     968             : /* same, public function */
     969             : GEN
     970           0 : RgM_diagonal(GEN m)
     971             : {
     972           0 :   long i, lx = lg(m);
     973           0 :   GEN y = cgetg(lx,t_VEC);
     974           0 :   for (i=1; i<lx; i++) gel(y,i) = gcopy(gcoeff(m,i,i));
     975           0 :   return y;
     976             : }
     977             : 
     978             : 

Generated by: LCOV version 1.11