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 - alglin3.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.16.2 lcov report (development 29419-8afb0ed749) Lines: 550 595 92.4 %
Date: 2024-07-02 09:03:41 Functions: 47 53 88.7 %
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; 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             : /********************************************************************/
      16             : /**                                                                **/
      17             : /**                         LINEAR ALGEBRA                         **/
      18             : /**                          (third part)                          **/
      19             : /**                                                                **/
      20             : /********************************************************************/
      21             : #include "pari.h"
      22             : #include "paripriv.h"
      23             : 
      24             : #define DEBUGLEVEL DEBUGLEVEL_mat
      25             : 
      26             : /*******************************************************************/
      27             : /*                                                                 */
      28             : /*                               SUM                               */
      29             : /*                                                                 */
      30             : /*******************************************************************/
      31             : 
      32             : GEN
      33      184543 : vecsum(GEN v)
      34             : {
      35      184543 :   pari_sp av = avma;
      36             :   long i, l;
      37             :   GEN p;
      38      184543 :   if (!is_vec_t(typ(v)))
      39           7 :     pari_err_TYPE("vecsum", v);
      40      184536 :   l = lg(v);
      41      184536 :   if (l == 1) return gen_0;
      42      184529 :   p = gel(v,1);
      43      184529 :   if (l == 2) return gcopy(p);
      44      248285 :   for (i=2; i<l; i++)
      45             :   {
      46      167865 :     p = gadd(p, gel(v,i));
      47      167865 :     if (gc_needed(av, 2))
      48             :     {
      49           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sum");
      50           0 :       p = gerepileupto(av, p);
      51             :     }
      52             :   }
      53       80420 :   return gerepileupto(av, p);
      54             : }
      55             : 
      56             : /*******************************************************************/
      57             : /*                                                                 */
      58             : /*                         TRANSPOSE                               */
      59             : /*                                                                 */
      60             : /*******************************************************************/
      61             : /* A[x0,]~ */
      62             : static GEN
      63    25477279 : row_transpose(GEN A, long x0)
      64             : {
      65    25477279 :   long i, lB = lg(A);
      66    25477279 :   GEN B  = cgetg(lB, t_COL);
      67   199991195 :   for (i=1; i<lB; i++) gel(B, i) = gcoeff(A, x0, i);
      68    25477293 :   return B;
      69             : }
      70             : static GEN
      71       18700 : row_transposecopy(GEN A, long x0)
      72             : {
      73       18700 :   long i, lB = lg(A);
      74       18700 :   GEN B  = cgetg(lB, t_COL);
      75      162385 :   for (i=1; i<lB; i++) gel(B, i) = gcopy(gcoeff(A, x0, i));
      76       18700 :   return B;
      77             : }
      78             : 
      79             : /* No copy*/
      80             : GEN
      81     6779215 : shallowtrans(GEN x)
      82             : {
      83             :   long i, dx, lx;
      84             :   GEN y;
      85     6779215 :   switch(typ(x))
      86             :   {
      87        1162 :     case t_VEC: y = leafcopy(x); settyp(y,t_COL); break;
      88         896 :     case t_COL: y = leafcopy(x); settyp(y,t_VEC); break;
      89     6777160 :     case t_MAT:
      90     6777160 :       lx = lg(x); if (lx==1) return cgetg(1,t_MAT);
      91     6775865 :       dx = lgcols(x); y = cgetg(dx,t_MAT);
      92    32253573 :       for (i = 1; i < dx; i++) gel(y,i) = row_transpose(x,i);
      93     6776457 :       break;
      94           0 :     default: pari_err_TYPE("shallowtrans",x);
      95             :       return NULL;/*LCOV_EXCL_LINE*/
      96             :   }
      97     6778515 :   return y;
      98             : }
      99             : 
     100             : GEN
     101       43617 : gtrans(GEN x)
     102             : {
     103             :   long i, dx, lx;
     104             :   GEN y;
     105       43617 :   switch(typ(x))
     106             :   {
     107       36771 :     case t_VEC: y = gcopy(x); settyp(y,t_COL); break;
     108        4879 :     case t_COL: y = gcopy(x); settyp(y,t_VEC); break;
     109        1960 :     case t_MAT:
     110        1960 :       lx = lg(x); if (lx==1) return cgetg(1,t_MAT);
     111        1953 :       dx = lgcols(x); y = cgetg(dx,t_MAT);
     112       20653 :       for (i = 1; i < dx; i++) gel(y,i) = row_transposecopy(x,i);
     113        1953 :       break;
     114           7 :     default: pari_err_TYPE("gtrans",x);
     115             :       return NULL;/*LCOV_EXCL_LINE*/
     116             :   }
     117       43603 :   return y;
     118             : }
     119             : 
     120             : /*******************************************************************/
     121             : /*                                                                 */
     122             : /*                           EXTRACTION                            */
     123             : /*                                                                 */
     124             : /*******************************************************************/
     125             : 
     126             : static long
     127         182 : str_to_long(char *s, char **pt)
     128             : {
     129         182 :   long a = atol(s);
     130         182 :   while (isspace((unsigned char)*s)) s++;
     131         182 :   if (*s == '-' || *s == '+') s++;
     132         385 :   while (isdigit((unsigned char)*s) || isspace((unsigned char)*s)) s++;
     133         182 :   *pt = s; return a;
     134             : }
     135             : 
     136             : static int
     137         112 : get_range(char *s, long *a, long *b, long *cmpl, long lx)
     138             : {
     139         112 :   long max = lx - 1;
     140             : 
     141         112 :   *a = 1; *b = max;
     142         112 :   if (*s == '^') { *cmpl = 1; s++; } else *cmpl = 0;
     143         112 :   if (!*s) return 0;
     144         112 :   if (*s != '.')
     145             :   {
     146         105 :     *a = str_to_long(s, &s);
     147         105 :     if (*a < 0) *a += lx;
     148         105 :     if (*a<1 || *a>max) return 0;
     149             :   }
     150         112 :   if (*s == '.')
     151             :   {
     152         105 :     s++; if (*s != '.') return 0;
     153         105 :     do s++; while (isspace((unsigned char)*s));
     154         105 :     if (*s)
     155             :     {
     156          77 :       *b = str_to_long(s, &s);
     157          77 :       if (*b < 0) *b += lx;
     158          77 :       if (*b<1 || *b>max || *s) return 0;
     159             :     }
     160          98 :     return 1;
     161             :   }
     162           7 :   if (*s) return 0;
     163           7 :   *b = *a; return 1;
     164             : }
     165             : 
     166             : static int
     167          35 : extract_selector_ok(long lx, GEN L)
     168             : {
     169             :   long i, l;
     170          35 :   switch (typ(L))
     171             :   {
     172           7 :     case t_INT: {
     173             :       long maxj;
     174           7 :       if (!signe(L)) return 1;
     175           7 :       l = lgefint(L)-1;
     176           7 :       maxj = BITS_IN_LONG - bfffo(*int_MSW(L));
     177           7 :       return ((l-2) * BITS_IN_LONG + maxj < lx);
     178             :     }
     179           7 :     case t_STR: {
     180             :       long first, last, cmpl;
     181           7 :       return get_range(GSTR(L), &first, &last, &cmpl, lx);
     182             :     }
     183          14 :     case t_VEC: case t_COL:
     184          14 :       l = lg(L);
     185          28 :       for (i=1; i<l; i++)
     186             :       {
     187          21 :         long j = itos(gel(L,i));
     188          21 :         if (j>=lx || j<=0) return 0;
     189             :       }
     190           7 :       return 1;
     191           7 :     case t_VECSMALL:
     192           7 :       l = lg(L);
     193          21 :       for (i=1; i<l; i++)
     194             :       {
     195          14 :         long j = L[i];
     196          14 :         if (j>=lx || j<=0) return 0;
     197             :       }
     198           7 :       return 1;
     199             :   }
     200           0 :   return 0;
     201             : }
     202             : 
     203             : GEN
     204       10297 : shallowmatextract(GEN x, GEN l1, GEN l2)
     205             : {
     206       10297 :   long i, j, n1 = lg(l1), n2 = lg(l2);
     207       10297 :   GEN M = cgetg(n2, t_MAT);
     208       65261 :   for(i=1; i < n2; i++)
     209             :   {
     210       54964 :     long ii = l2[i];
     211       54964 :     GEN C = cgetg(n1, t_COL);
     212      957404 :     for (j=1; j < n1; j++)
     213             :     {
     214      902440 :       long jj = l1[j];
     215      902440 :       gel(C, j) = gmael(x, ii, jj);
     216             :     }
     217       54964 :     gel(M, i) = C;
     218             :   }
     219       10297 :   return M;
     220             : }
     221             : 
     222             : GEN
     223       39101 : shallowextract(GEN x, GEN L)
     224             : {
     225       39101 :   long i,j, tl = typ(L), tx = typ(x), lx = lg(x);
     226             :   GEN y;
     227             : 
     228       39101 :   switch(tx)
     229             :   {
     230       39094 :     case t_VEC:
     231             :     case t_COL:
     232             :     case t_MAT:
     233       39094 :     case t_VECSMALL: break;
     234           7 :     default: pari_err_TYPE("extract",x);
     235             : 
     236             :   }
     237       39094 :   if (tl==t_INT)
     238             :   { /* extract components of x as per the bits of mask L */
     239             :     long k, l, ix, iy, maxj;
     240             :     GEN Ld;
     241        3318 :     if (!signe(L)) return cgetg(1,tx);
     242        3311 :     y = new_chunk(lx);
     243        3311 :     l = lgefint(L)-1; ix = iy = 1;
     244        3311 :     maxj = BITS_IN_LONG - bfffo(*int_MSW(L));
     245        3311 :     if ((l-2) * BITS_IN_LONG + maxj >= lx)
     246           7 :       pari_err_TYPE("vecextract [mask too large]", L);
     247        3679 :     for (k = 2, Ld = int_LSW(L); k < l; k++, Ld = int_nextW(Ld))
     248             :     {
     249         375 :       ulong B = *Ld;
     250       20439 :       for (j = 0; j < BITS_IN_LONG; j++, B >>= 1, ix++)
     251       20064 :         if (B & 1) y[iy++] = x[ix];
     252             :     }
     253             :     { /* k = l */
     254        3304 :       ulong B = *Ld;
     255       28698 :       for (j = 0; j < maxj; j++, B >>= 1, ix++)
     256       25394 :         if (B & 1) y[iy++] = x[ix];
     257             :     }
     258        3304 :     y[0] = evaltyp(tx) | evallg(iy);
     259        3304 :     return y;
     260             :   }
     261       35776 :   if (tl==t_STR)
     262             :   {
     263         105 :     char *s = GSTR(L);
     264             :     long first, last, cmpl, d;
     265         105 :     if (! get_range(s, &first, &last, &cmpl, lx))
     266           7 :       pari_err_TYPE("vecextract [incorrect range]", L);
     267          98 :     if (lx == 1) return cgetg(1,tx);
     268          98 :     d = last - first;
     269          98 :     if (cmpl)
     270             :     {
     271          21 :       if (d >= 0)
     272             :       {
     273          14 :         y = cgetg(lx - (1+d),tx);
     274         469 :         for (j=1; j<first; j++) gel(y,j) = gel(x,j);
     275         266 :         for (i=last+1; i<lx; i++,j++) gel(y,j) = gel(x,i);
     276             :       }
     277             :       else
     278             :       {
     279           7 :         y = cgetg(lx - (1-d),tx);
     280          14 :         for (j=1,i=lx-1; i>first; i--,j++) gel(y,j) = gel(x,i);
     281          14 :         for (i=last-1; i>0; i--,j++) gel(y,j) = gel(x,i);
     282             :       }
     283             :     }
     284             :     else
     285             :     {
     286          77 :       if (d >= 0)
     287             :       {
     288          35 :         y = cgetg(d+2,tx);
     289         112 :         for (i=first,j=1; i<=last; i++,j++) gel(y,j) = gel(x,i);
     290             :       }
     291             :       else
     292             :       {
     293          42 :         y = cgetg(2-d,tx);
     294         203 :         for (i=first,j=1; i>=last; i--,j++) gel(y,j) = gel(x,i);
     295             :       }
     296             :     }
     297          98 :     return y;
     298             :   }
     299             : 
     300       35671 :   if (is_vec_t(tl))
     301             :   {
     302          77 :     long ll=lg(L); y=cgetg(ll,tx);
     303         196 :     for (i=1; i<ll; i++)
     304             :     {
     305         133 :       j = itos(gel(L,i));
     306         133 :       if (j<=0) pari_err_COMPONENT("vecextract","<=",gen_0,stoi(j));
     307         126 :       if (j>=lx) pari_err_COMPONENT("vecextract",">=",stoi(lx),stoi(j));
     308         119 :       gel(y,i) = gel(x,j);
     309             :     }
     310          63 :     return y;
     311             :   }
     312       35594 :   if (tl == t_VECSMALL)
     313             :   {
     314       35587 :     long ll=lg(L); y=cgetg(ll,tx);
     315      159061 :     for (i=1; i<ll; i++)
     316             :     {
     317      123474 :       j = L[i];
     318      123474 :       if (j<=0) pari_err_COMPONENT("vecextract","<=",gen_0,stoi(j));
     319      123474 :       if (j>=lx) pari_err_COMPONENT("vecextract",">=",stoi(lx),stoi(j));
     320      123474 :       gel(y,i) = gel(x,j);
     321             :     }
     322       35587 :     return y;
     323             :   }
     324           7 :   pari_err_TYPE("vecextract [mask]", L);
     325             :   return NULL; /* LCOV_EXCL_LINE */
     326             : }
     327             : 
     328             : /* does the component selector l select 0 component ? */
     329             : static int
     330          77 : select_0(GEN l)
     331             : {
     332          77 :   switch(typ(l))
     333             :   {
     334          14 :     case t_INT:
     335          14 :       return (!signe(l));
     336          42 :     case t_VEC: case t_COL: case t_VECSMALL:
     337          42 :       return (lg(l) == 1);
     338             :   }
     339          21 :   return 0;
     340             : }
     341             : 
     342             : GEN
     343       29111 : extract0(GEN x, GEN l1, GEN l2)
     344             : {
     345       29111 :   pari_sp av = avma, av2;
     346             :   GEN y;
     347       29111 :   if (! l2)
     348             :   {
     349       29034 :     y = shallowextract(x, l1);
     350       28992 :     if (lg(y) == 1 || typ(y) == t_VECSMALL) return y;
     351       28985 :     av2 = avma;
     352       28985 :     y = gcopy(y);
     353             :   }
     354             :   else
     355             :   {
     356          77 :     if (typ(x) != t_MAT) pari_err_TYPE("extract",x);
     357          77 :     y = shallowextract(x,l2);
     358          77 :     if (select_0(l1)) { set_avma(av); return zeromat(0, lg(y)-1); }
     359          63 :     if (lg(y) == 1 && lg(x) > 1)
     360             :     {
     361          35 :       if (!extract_selector_ok(lgcols(x), l1))
     362           7 :         pari_err_TYPE("vecextract [incorrect mask]", l1);
     363          28 :       set_avma(av); return cgetg(1, t_MAT);
     364             :     }
     365          28 :     y = shallowextract(shallowtrans(y), l1);
     366          28 :     av2 = avma;
     367          28 :     y = gtrans(y);
     368             :   }
     369       29013 :   stackdummy(av, av2);
     370       29013 :   return y;
     371             : }
     372             : 
     373             : static long
     374        2016 : vecslice_parse_arg(long lA, long *y1, long *y2, long *skip)
     375             : {
     376        2016 :   *skip=0;
     377        2016 :   if (*y1==LONG_MAX)
     378             :   {
     379         252 :     if (*y2!=LONG_MAX)
     380             :     {
     381         140 :       if (*y2<0) *y2 += lA;
     382         140 :       if (*y2<0 || *y2==LONG_MAX || *y2>=lA)
     383           0 :         pari_err_DIM("_[..]");
     384         140 :       *skip=*y2;
     385             :     }
     386         252 :     *y1 = 1; *y2 = lA-1;
     387             :   }
     388        1764 :   else if (*y2==LONG_MAX) *y2 = *y1;
     389        2016 :   if (*y1<=0) *y1 += lA;
     390        2016 :   if (*y2<0) *y2 += lA;
     391        2016 :   if (*y1<=0 || *y1>*y2+1 || *y2>=lA) pari_err_DIM("_[..]");
     392        2002 :   return *y2 - *y1 + 2 - !!*skip;
     393             : }
     394             : 
     395             : static GEN
     396        2653 : vecslice_i(GEN A, long t, long lB, long y1, long skip)
     397             : {
     398        2653 :   GEN B = cgetg(lB, t);
     399             :   long i;
     400       27034 :   for (i=1; i<lB; i++, y1++)
     401             :   {
     402       24381 :     if (y1 == skip) { i--; continue; }
     403       24241 :     gel(B,i) = gcopy(gel(A,y1));
     404             :   }
     405        2653 :   return B;
     406             : }
     407             : 
     408             : static GEN
     409          14 : rowslice_i(GEN A, long lB, long x1, long y1, long skip)
     410             : {
     411          14 :   GEN B = cgetg(lB, t_VEC);
     412             :   long i;
     413          77 :   for (i=1; i<lB; i++, y1++)
     414             :   {
     415          63 :     if (y1 == skip) { i--; continue; }
     416          56 :     gel(B,i) = gcopy(gcoeff(A,x1,y1));
     417             :   }
     418          14 :   return B;
     419             : }
     420             : 
     421             : static GEN
     422           0 : rowsmallslice_i(GEN A, long lB, long x1, long y1, long skip)
     423             : {
     424           0 :   GEN B = cgetg(lB, t_VECSMALL);
     425             :   long i;
     426           0 :   for (i=1; i<lB; i++, y1++)
     427             :   {
     428           0 :     if (y1 == skip) { i--; continue; }
     429           0 :     B[i] = coeff(A,x1,y1);
     430             :   }
     431           0 :   return B;
     432             : }
     433             : 
     434             : static GEN
     435          28 : vecsmallslice_i(GEN A, long t, long lB, long y1, long skip)
     436             : {
     437          28 :   GEN B = cgetg(lB, t);
     438             :   long i;
     439         126 :   for (i=1; i<lB; i++, y1++)
     440             :   {
     441          98 :     if (y1 == skip) { i--; continue; }
     442          91 :     B[i] = A[y1];
     443             :   }
     444          28 :   return B;
     445             : }
     446             : GEN
     447        1617 : vecslice0(GEN A, long y1, long y2)
     448             : {
     449        1617 :   long skip, lB, t = typ(A);
     450        1617 :   switch(t)
     451             :   {
     452        1519 :     case t_VEC: case t_COL:
     453        1519 :       lB = vecslice_parse_arg(lg(A), &y1, &y2, &skip);
     454        1505 :       return vecslice_i(A, t,lB,y1,skip);
     455          28 :     case t_VECSMALL:
     456          28 :       lB = vecslice_parse_arg(lg(A), &y1, &y2, &skip);
     457          28 :       return vecsmallslice_i(A, t,lB,y1,skip);
     458          63 :     case t_LIST:
     459          63 :       if (list_typ(A) == t_LIST_RAW)
     460             :       {
     461          63 :         GEN y, z = list_data(A);
     462          63 :         long l = z? lg(z): 1;
     463          63 :         lB = vecslice_parse_arg(l, &y1, &y2, &skip);
     464          63 :         y = mklist(); if (!z) return y;
     465          63 :         list_data(y) = vecslice_i(z, t_VEC,lB,y1,skip);
     466          63 :         return y;
     467             :       }
     468             :     default:
     469           7 :       pari_err_TYPE("_[_.._]",A);
     470             :       return NULL;/*LCOV_EXCL_LINE*/
     471             :   }
     472             : }
     473             : 
     474             : GEN
     475         210 : matslice0(GEN A, long x1, long x2, long y1, long y2)
     476             : {
     477             :   GEN B;
     478         210 :   long i, lB, lA = lg(A), rA, t, skip, rskip, rlB;
     479         210 :   long is_col = y1!=LONG_MAX && y2==LONG_MAX;
     480         210 :   long is_row = x1!=LONG_MAX && x2==LONG_MAX;
     481             :   GEN (*slice)(GEN A, long t, long lB, long y1, long skip);
     482         210 :   if (typ(A)!=t_MAT) pari_err_TYPE("_[_.._,_.._]",A);
     483         210 :   lB = vecslice_parse_arg(lA, &y1, &y2, &skip);
     484         210 :   if (is_col) return vecslice0(gel(A, y1), x1, x2);
     485         196 :   rA = lg(A)==1 ? 1: lgcols(A);
     486         196 :   rlB = vecslice_parse_arg(rA, &x1, &x2, &rskip);
     487         196 :   t = lg(A)==1 ? t_COL: typ(gel(A,1));
     488         196 :   if (is_row) return t == t_COL ? rowslice_i(A, lB, x1, y1, skip):
     489           0 :                                   rowsmallslice_i(A, lB, x1, y1, skip);
     490         182 :   slice = t == t_COL? &vecslice_i: &vecsmallslice_i;
     491             : 
     492         182 :   B = cgetg(lB, t_MAT);
     493        1281 :   for (i=1; i<lB; i++, y1++)
     494             :   {
     495        1099 :     if (y1 == skip) { i--; continue; }
     496        1085 :     gel(B,i) = slice(gel(A,y1),t,rlB, x1, rskip);
     497             :   }
     498         182 :   return B;
     499             : }
     500             : 
     501             : GEN
     502       10660 : vecrange(GEN a, GEN b)
     503             : {
     504             :   GEN y;
     505             :   long i, l;
     506       10660 :   if (typ(a)!=t_INT) pari_err_TYPE("[_.._]",a);
     507       10653 :   if (typ(b)!=t_INT) pari_err_TYPE("[_.._]",b);
     508       10646 :   if (cmpii(a,b)>0) return cgetg(1,t_VEC);
     509       10639 :   l = itos(subii(b,a))+1;
     510       10639 :   a = setloop(a);
     511       10639 :   y = cgetg(l+1, t_VEC);
     512    26384689 :   for (i=1; i<=l; a = incloop(a), i++) gel(y,i) = icopy(a);
     513       10639 :   return y;
     514             : }
     515             : 
     516             : GEN
     517           0 : vecrangess(long a, long b)
     518             : {
     519             :   GEN y;
     520             :   long i, l;
     521           0 :   if (a>b) return cgetg(1,t_VEC);
     522           0 :   l = b-a+1;
     523           0 :   y = cgetg(l+1, t_VEC);
     524           0 :   for (i=1; i<=l; a++, i++) gel(y,i) = stoi(a);
     525           0 :   return y;
     526             : }
     527             : 
     528             : GEN
     529         110 : genindexselect(void *E, long (*f)(void* E, GEN x), GEN A)
     530             : {
     531             :   long l, i, lv;
     532             :   GEN v, z;
     533             :   pari_sp av;
     534         110 :   switch(typ(A))
     535             :   {
     536          14 :     case t_LIST:
     537          14 :       z = list_data(A);
     538          14 :       l = z? lg(z): 1;
     539          14 :       if (list_typ(A)==t_LIST_MAP)
     540             :       {
     541           7 :         av = avma;
     542           7 :         return gerepilecopy(av, mapselect_shallow(E, f, A));
     543             :       }
     544           7 :       break;
     545          89 :     case t_VEC: case t_COL: case t_MAT:
     546          89 :       l = lg(A);
     547          89 :       z = A;
     548          89 :       break;
     549           7 :     default:
     550           7 :       pari_err_TYPE("select",A);
     551             :       return NULL;/*LCOV_EXCL_LINE*/
     552             :   }
     553          96 :   v = cgetg(l, t_VECSMALL);
     554          96 :   av = avma;
     555          96 :   clone_lock(A);
     556       12878 :   for (i = lv = 1; i < l; i++) {
     557       12782 :     if (f(E, gel(z,i))) v[lv++] = i;
     558       12782 :     set_avma(av);
     559             :   }
     560          96 :   clone_unlock_deep(A); fixlg(v, lv); return v;
     561             : }
     562             : static GEN
     563         101 : extract_copy(GEN A, GEN v)
     564             : {
     565         101 :   long i, l = lg(v);
     566         101 :   GEN B = cgetg(l, typ(A));
     567        4609 :   for (i = 1; i < l; i++) gel(B,i) = gcopy(gel(A,v[i]));
     568         101 :   return B;
     569             : }
     570             : /* as genselect, but treat A [ t_VEC,t_COL, or t_MAT] as a t_VEC */
     571             : GEN
     572           0 : vecselect(void *E, long (*f)(void* E, GEN x), GEN A)
     573             : {
     574             :   GEN v;
     575           0 :   clone_lock(A);
     576           0 :   v = genindexselect(E, f, A);
     577           0 :   A = extract_copy(A, v); settyp(A, t_VEC);
     578           0 :   clone_unlock_deep(A); return A;
     579             : }
     580             : GEN
     581         104 : genselect(void *E, long (*f)(void* E, GEN x), GEN A)
     582             : {
     583         104 :   pari_sp av  = avma;
     584             :   GEN y, z, v;/* v left on stack for efficiency */
     585         104 :   clone_lock(A);
     586         104 :   switch(typ(A))
     587             :   {
     588          35 :     case t_LIST:
     589          35 :       z = list_data(A);
     590          35 :       if (!z) y = mklist();
     591             :       else
     592             :       {
     593          35 :         if (list_typ(A)==t_LIST_MAP)
     594             :         {
     595          14 :           long i, l = z? lg(z): 1, lv=1;
     596          14 :           GEN v1 = cgetg(l, t_COL);
     597          14 :           GEN v2 = cgetg(l, t_COL);
     598          56 :           for (i = lv = 1; i < l; i++) {
     599          42 :             if (f(E, gmael3(z,i,1,2)))
     600             :            {
     601          21 :              gel(v1,lv) = gmael3(z,i,1,1);
     602          21 :              gel(v2,lv) = gmael3(z,i,1,2);
     603          21 :              lv++;
     604             :            }
     605             :           }
     606          14 :           fixlg(v1, lv); fixlg(v2, lv); y = gtomap(mkmat2(v1,v2));
     607             :         }
     608             :         else
     609             :         {
     610             :           GEN B;
     611          21 :           y = cgetg(3, t_LIST);
     612          21 :           v = genindexselect(E, f, z);
     613          21 :           B = extract_copy(z, v);
     614          21 :           y[1] = lg(B)-1;
     615          21 :           list_data(y) = B;
     616             :         }
     617             :       }
     618          35 :       break;
     619          62 :     case t_VEC: case t_COL: case t_MAT:
     620          62 :       v = genindexselect(E, f, A);
     621          62 :       y = extract_copy(A, v);
     622          62 :       break;
     623           7 :     default:
     624           7 :       pari_err_TYPE("select",A);
     625             :       return NULL;/*LCOV_EXCL_LINE*/
     626             :   }
     627          97 :   clone_unlock_deep(A); return gerepileupto(av, y);
     628             : }
     629             : 
     630             : static void
     631       54731 : check_callgen1(GEN f, const char *s)
     632             : {
     633       54731 :   if (typ(f) != t_CLOSURE || closure_is_variadic(f)  || closure_arity(f) < 1)
     634           0 :     pari_err_TYPE(s, f);
     635       54731 : }
     636             : 
     637             : GEN
     638         131 : select0(GEN f, GEN x, long flag)
     639             : {
     640         131 :   check_callgen1(f, "select");
     641         131 :   switch(flag)
     642             :   {
     643         104 :     case 0: return genselect((void *) f, gp_callbool, x);
     644          27 :     case 1: return genindexselect((void *) f, gp_callbool, x);
     645           0 :     default: pari_err_FLAG("select");
     646             :              return NULL;/*LCOV_EXCL_LINE*/
     647             :   }
     648             : }
     649             : 
     650             : GEN
     651          30 : parselect(GEN C, GEN D, long flag)
     652             : {
     653             :   pari_sp av, av2;
     654          30 :   long lv, l = lg(D), i, pending = 0, workid;
     655             :   GEN V, done;
     656             :   struct pari_mt pt;
     657          30 :   check_callgen1(C, "parselect");
     658          30 :   if (!is_vec_t(typ(D))) pari_err_TYPE("parselect",D);
     659          30 :   V = cgetg(l, t_VECSMALL); av = avma;
     660          30 :   mt_queue_start_lim(&pt, C, l-1);
     661          30 :   av2 = avma;
     662       30204 :   for (i=1; i<l || pending; i++)
     663             :   {
     664       30174 :     mt_queue_submit(&pt, i, i<l? mkvec(gel(D,i)): NULL);
     665       30174 :     done = mt_queue_get(&pt, &workid, &pending);
     666       30174 :     if (done) V[workid] = !gequal0(done);
     667       30174 :     set_avma(av2);
     668             :   }
     669          30 :   mt_queue_end(&pt);
     670          30 :   set_avma(av);
     671       30054 :   for (lv=1, i=1; i<l; i++)
     672       30024 :     if (V[i]) V[lv++]=i;
     673          30 :   fixlg(V, lv);
     674          30 :   return flag? V: extract_copy(D, V);
     675             : }
     676             : 
     677             : GEN
     678           0 : veccatapply(void *E, GEN (*f)(void* E, GEN x), GEN x)
     679             : {
     680           0 :   pari_sp av = avma;
     681           0 :   GEN v = vecapply(E, f, x);
     682           0 :   return lg(v) == 1? v: gerepilecopy(av, shallowconcat1(v));
     683             : }
     684             : 
     685             : static GEN
     686          14 : vecapply2(void *E, GEN (*f)(void* E, GEN x), GEN x)
     687             : {
     688             :   long i, lx;
     689          14 :   GEN y = cgetg_copy(x, &lx); y[1] = x[1];
     690          56 :   for (i=2; i<lx; i++) gel(y,i) = f(E, gel(x,i));
     691          14 :   return y;
     692             : }
     693             : static GEN
     694      168232 : vecapply1(void *E, GEN (*f)(void* E, GEN x), GEN x)
     695             : {
     696             :   long i, lx;
     697      168232 :   GEN y = cgetg_copy(x, &lx);
     698      864668 :   for (i=1; i<lx; i++) gel(y,i) = f(E, gel(x,i));
     699      168232 :   return y;
     700             : }
     701             : static GEN
     702          21 : mapapply1(void *E, GEN (*f)(void* E, GEN x), GEN x)
     703             : {
     704             :   long i, lx;
     705          21 :   GEN y = cgetg_copy(x, &lx);
     706          84 :   for (i=1; i<lx; i++)
     707             :   {
     708          63 :     GEN xi = gel(x, i);
     709          63 :     gel(y,i) = mkvec2(mkvec2(gcopy(gmael(xi, 1, 1)), f(E,gmael(xi, 1, 2))),
     710          63 :                              gcopy(gel(xi, 2)));
     711             :   }
     712          21 :   return y;
     713             : }
     714             : /* as genapply, but treat A [ t_VEC,t_COL, or t_MAT] as a t_VEC */
     715             : GEN
     716      115406 : vecapply(void *E, GEN (*f)(void* E, GEN x), GEN x)
     717             : {
     718             :   GEN y;
     719      115406 :   clone_lock(x); y = vecapply1(E,f,x);
     720      115406 :   clone_unlock_deep(x); settyp(y, t_VEC); return y;
     721             : }
     722             : GEN
     723       52826 : genapply(void *E, GEN (*f)(void* E, GEN x), GEN x)
     724             : {
     725       52826 :   long i, lx, tx = typ(x);
     726             :   GEN y, z;
     727       52826 :   if (is_scalar_t(tx)) return f(E, x);
     728       52826 :   clone_lock(x);
     729       52826 :   switch(tx) {
     730           7 :     case t_POL: y = normalizepol(vecapply2(E,f,x)); break;
     731           7 :     case t_SER:
     732           7 :       y = ser_isexactzero(x)? gcopy(x): normalizeser(vecapply2(E,f,x));
     733           7 :       break;
     734          42 :     case t_LIST:
     735             :       {
     736          42 :         long t = list_typ(x);
     737          42 :         z = list_data(x);
     738          42 :         if (!z)
     739           7 :           y = mklist_typ(t);
     740             :         else
     741             :         {
     742          35 :           y = cgetg(3, t_LIST);
     743          35 :           y[1] = evaltyp(t)|_evallg(lg(z)-1);
     744             :           switch(t)
     745             :           {
     746          14 :           case t_LIST_RAW:
     747          14 :             list_data(y) = vecapply1(E,f,z);
     748          14 :             break;
     749          21 :           case t_LIST_MAP:
     750          21 :             list_data(y) = mapapply1(E,f,z);
     751          21 :             break;
     752             :           }
     753             :         }
     754             :       }
     755          42 :       break;
     756          42 :     case t_MAT:
     757          42 :       y = cgetg_copy(x, &lx);
     758         126 :       for (i = 1; i < lx; i++) gel(y,i) = vecapply1(E,f,gel(x,i));
     759          42 :       break;
     760             : 
     761       52728 :     case t_VEC: case t_COL: y = vecapply1(E,f,x); break;
     762           0 :     default:
     763           0 :       pari_err_TYPE("apply",x);
     764             :       return NULL;/*LCOV_EXCL_LINE*/
     765             :   }
     766       52826 :   clone_unlock_deep(x); return y;
     767             : }
     768             : 
     769             : GEN
     770       52826 : apply0(GEN f, GEN x)
     771             : {
     772       52826 :   check_callgen1(f, "apply");
     773       52826 :   return genapply((void *) f, gp_call, x);
     774             : }
     775             : 
     776             : GEN
     777         448 : vecselapply(void *Epred, long (*pred)(void* E, GEN x), void *Efun,
     778             :                          GEN (*fun)(void* E, GEN x), GEN A)
     779             : {
     780             :   GEN y;
     781         448 :   long i, l = lg(A), nb=1;
     782         448 :   clone_lock(A); y = cgetg(l, t_VEC);
     783    26165475 :   for (i=1; i<l; i++)
     784    26165027 :     if (pred(Epred, gel(A,i))) gel(y,nb++) = fun(Efun, gel(A,i));
     785         448 :   fixlg(y,nb); clone_unlock_deep(A); return y;
     786             : }
     787             : 
     788             : GEN
     789           0 : veccatselapply(void *Epred, long (*pred)(void* E, GEN x), void *Efun,
     790             :                             GEN (*fun)(void* E, GEN x), GEN A)
     791             : {
     792           0 :   pari_sp av = avma;
     793           0 :   GEN v = vecselapply(Epred, pred, Efun, fun, A);
     794           0 :   return lg(v) == 1? v: gerepilecopy(av, shallowconcat1(v));
     795             : }
     796             : 
     797             : /* suitable for gerepileupto */
     798             : GEN
     799          44 : parapply_slice_worker(GEN D, GEN worker)
     800             : {
     801             :   long l, i;
     802          44 :   GEN w = cgetg_copy(D, &l);
     803       17522 :   for (i = 1; i < l; i++) gel(w,i) = closure_callgen1(worker, gel(D,i));
     804          43 :   return w;
     805             : }
     806             : 
     807             : /* B <- {A[i] : i = r (mod m)}, 1 <= r <= m */
     808             : static void
     809       73185 : arithprogset(GEN B, GEN A, long r, long m)
     810             : {
     811       73185 :   long i, k, l = lg(A);
     812      170775 :   for (k = 1, i = r; i < l; i += m, k++) gel(B, k) = gel(A, i);
     813       73185 :   setlg(B, k);
     814       73185 : }
     815             : GEN
     816        7246 : gen_parapply_slice(GEN worker, GEN D, long mmin)
     817             : {
     818        7246 :   long l, r, n = lg(D)-1, m = minss(mmin, n), pending = 0;
     819        7246 :   GEN L = cgetg(n / m + 2, t_VEC), va = mkvec(L), V = cgetg_copy(D, &l);
     820             :   struct pari_mt pt;
     821        7246 :   mt_queue_start_lim(&pt, worker, m);
     822      146370 :   for (r = 1; r <= m || pending; r++)
     823             :   {
     824             :     long workid;
     825             :     GEN done;
     826      139124 :     if (r <= m) arithprogset(L, D, r, m);
     827      139124 :     mt_queue_submit(&pt, r, r <= m? va: NULL);
     828      139124 :     done = mt_queue_get(&pt, &workid, &pending);
     829      139124 :     if (done)
     830             :     {
     831       73185 :       long j, k, J = lg(done)-1;
     832      170775 :       for (j = 1, k = workid; j <= J; j++, k +=m) gel(V, k) = gel(done, j);
     833             :     }
     834             :   }
     835        7246 :   mt_queue_end(&pt); return V;
     836             : }
     837             : 
     838             : GEN
     839      374376 : gen_parapply_percent(GEN worker, GEN D, long percent)
     840             : {
     841      374376 :   long l = lg(D), i, pending = 0, cnt = 0, lper = -1, lcnt = 0;
     842             :   GEN W, V;
     843             :   struct pari_mt pt;
     844             : 
     845      374376 :   if (l == 1) return cgetg(1, typ(D));
     846      343954 :   W = cgetg(2, t_VEC);
     847      343952 :   V = cgetg(l, typ(D));
     848      343952 :   mt_queue_start_lim(&pt, worker, l-1);
     849     2751812 :   for (i = 1; i < l || pending; i++)
     850             :   {
     851             :     long workid;
     852             :     GEN done;
     853     2407864 :     if (i < l) gel(W,1) = gel(D,i);
     854     2407864 :     mt_queue_submit(&pt, i, i<l? W: NULL);
     855     2407812 :     done = mt_queue_get(&pt, &workid, &pending);
     856     2407802 :     if (done)
     857             :     {
     858     2345989 :       gel(V,workid) = done;
     859     2345989 :       if (percent && (++cnt)-lcnt>=percent)
     860             :       {
     861           0 :         long per = (long)(cnt*100./(l-1));
     862           0 :         lcnt = cnt;
     863           0 :         if (per > lper) { err_printf("%ld%% ",per); lper = per; }
     864             :       }
     865             :     }
     866             :   }
     867      343948 :   mt_queue_end(&pt); return V;
     868             : }
     869             : 
     870             : GEN
     871      371471 : gen_parapply(GEN worker, GEN D)
     872             : {
     873      371471 :   return gen_parapply_percent(worker, D, 0);
     874             : }
     875             : 
     876             : GEN
     877        1744 : parapply(GEN C, GEN D)
     878             : {
     879        1744 :   pari_sp av = avma;
     880        1744 :   check_callgen1(C, "parapply");
     881        1744 :   if (!is_vec_t(typ(D))) pari_err_TYPE("parapply",D);
     882        1744 :   return gerepileupto(av, gen_parapply(C, D));
     883             : }
     884             : 
     885             : GEN
     886          28 : genfold(void *E, GEN (*f)(void* E, GEN x, GEN y), GEN x)
     887             : {
     888          28 :   pari_sp av = avma;
     889             :   GEN z;
     890          28 :   long i, l = lg(x);
     891          28 :   if (!is_vec_t(typ(x))|| l==1  ) pari_err_TYPE("fold",x);
     892          28 :   clone_lock(x);
     893          28 :   z = gel(x,1);
     894         119 :   for (i=2; i<l; i++)
     895             :   {
     896          91 :     z = f(E,z,gel(x,i));
     897          91 :     if (gc_needed(av, 2))
     898             :     {
     899           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"fold");
     900           0 :       z = gerepilecopy(av, z);
     901             :     }
     902             :   }
     903          28 :   clone_unlock_deep(x);
     904          28 :   return gerepilecopy(av, z);
     905             : }
     906             : 
     907             : GEN
     908          28 : fold0(GEN f, GEN x)
     909             : {
     910          28 :   if (typ(f) != t_CLOSURE || closure_arity(f) < 2) pari_err_TYPE("fold",f);
     911          28 :   return genfold((void *) f, gp_call2, x);
     912             : }
     913             : /*******************************************************************/
     914             : /*                                                                 */
     915             : /*                     SCALAR-MATRIX OPERATIONS                    */
     916             : /*                                                                 */
     917             : /*******************************************************************/
     918             : GEN
     919      347663 : gtomat(GEN x)
     920             : {
     921             :   long lx, i;
     922             :   GEN y;
     923             : 
     924      347663 :   if (!x) return cgetg(1, t_MAT);
     925      347614 :   switch(typ(x))
     926             :   {
     927          28 :     case t_LIST:
     928          28 :       if (list_typ(x)==t_LIST_MAP)
     929          14 :         return maptomat(x);
     930          14 :       x = list_data(x);
     931          14 :       if (!x) return cgetg(1, t_MAT);
     932             :       /* fall through */
     933             :     case t_VEC: {
     934        2783 :       lx=lg(x); y=cgetg(lx,t_MAT);
     935        2786 :       if (lx == 1) break;
     936        2786 :       if (typ(gel(x,1)) == t_COL) {
     937        2415 :         long h = lgcols(x);
     938        5698 :         for (i=2; i<lx; i++) {
     939        3283 :           if (typ(gel(x,i)) != t_COL || lg(gel(x,i)) != h) break;
     940             :         }
     941        2415 :         if (i == lx) { /* matrix with h-1 rows */
     942        2415 :           y = cgetg(lx, t_MAT);
     943        8113 :           for (i=1 ; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
     944        2415 :           return y;
     945             :         }
     946             :       }
     947        1106 :       for (i=1; i<lx; i++) gel(y,i) = mkcolcopy(gel(x,i));
     948         371 :       break;
     949             :     }
     950      112209 :     case t_COL:
     951      112209 :       lx = lg(x);
     952      112209 :       if (lx == 1) return cgetg(1, t_MAT);
     953      112195 :       if (typ(gel(x,1)) == t_VEC) {
     954           7 :         long j, h = lg(gel(x,1));
     955          14 :         for (i=2; i<lx; i++) {
     956           7 :           if (typ(gel(x,i)) != t_VEC || lg(gel(x,i)) != h) break;
     957             :         }
     958           7 :         if (i == lx) { /* matrix with h cols */
     959           7 :           y = cgetg(h, t_MAT);
     960          28 :           for (j=1 ; j<h; j++) {
     961          21 :             gel(y,j) = cgetg(lx, t_COL);
     962          63 :             for (i=1; i<lx; i++) gcoeff(y,i,j) = gcopy(gmael(x,i,j));
     963             :           }
     964           7 :           return y;
     965             :         }
     966             :       }
     967      112188 :       y = mkmatcopy(x); break;
     968      231684 :     case t_MAT:
     969      231684 :       y = gcopy(x); break;
     970           7 :     case t_QFB: {
     971             :       GEN b;
     972           7 :       y = cgetg(3,t_MAT); b = gmul2n(gel(x,2),-1);
     973           7 :       gel(y,1) = mkcol2(icopy(gel(x,1)), b);
     974           7 :       gel(y,2) = mkcol2(b, icopy(gel(x,3)));
     975           7 :       break;
     976             :     }
     977         910 :     default:
     978         910 :       y = cgetg(2,t_MAT); gel(y,1) = mkcolcopy(x);
     979         910 :       break;
     980             :   }
     981      345162 :   return y;
     982             : }
     983             : 
     984             : /* create the diagonal matrix, whose diagonal is given by x */
     985             : GEN
     986     1334581 : diagonal(GEN x)
     987             : {
     988     1334581 :   long j, lx, tx = typ(x);
     989             :   GEN y;
     990             : 
     991     1334581 :   if (! is_matvec_t(tx)) return scalarmat(x,1);
     992     1334574 :   if (tx==t_MAT)
     993             :   {
     994          14 :     if (RgM_isdiagonal(x)) return gcopy(x);
     995           7 :     pari_err_TYPE("diagonal",x);
     996             :   }
     997     1334560 :   lx=lg(x); y=cgetg(lx,t_MAT);
     998     4761613 :   for (j=1; j<lx; j++)
     999             :   {
    1000     3427053 :     gel(y,j) = zerocol(lx-1);
    1001     3427053 :     gcoeff(y,j,j) = gcopy(gel(x,j));
    1002             :   }
    1003     1334560 :   return y;
    1004             : }
    1005             : /* same, assuming x is a t_VEC/t_COL. Not memory clean. */
    1006             : GEN
    1007      350455 : diagonal_shallow(GEN x)
    1008             : {
    1009      350455 :   long j, lx = lg(x);
    1010      350455 :   GEN y = cgetg(lx,t_MAT);
    1011             : 
    1012      962926 :   for (j=1; j<lx; j++)
    1013             :   {
    1014      612436 :     gel(y,j) = zerocol(lx-1);
    1015      612451 :     gcoeff(y,j,j) = gel(x,j);
    1016             :   }
    1017      350490 :   return y;
    1018             : }
    1019             : 
    1020             : GEN
    1021         448 : zv_diagonal(GEN x)
    1022             : {
    1023         448 :   long j, l = lg(x), n = l-1;
    1024         448 :   GEN y = cgetg(l,t_MAT);
    1025             : 
    1026        1428 :   for (j = 1; j < l; j++)
    1027             :   {
    1028         980 :     gel(y,j) = zero_Flv(n);
    1029         980 :     ucoeff(y,j,j) = uel(x,j);
    1030             :   }
    1031         448 :   return y;
    1032             : }
    1033             : 
    1034             : /* compute m*diagonal(d) */
    1035             : GEN
    1036          70 : matmuldiagonal(GEN m, GEN d)
    1037             : {
    1038             :   long j, lx;
    1039          70 :   GEN y = cgetg_copy(m, &lx);
    1040             : 
    1041          70 :   if (typ(m)!=t_MAT) pari_err_TYPE("matmuldiagonal",m);
    1042          70 :   if (! is_vec_t(typ(d))) pari_err_TYPE("matmuldiagonal",d);
    1043          70 :   if (lg(d) != lx) pari_err_OP("operation 'matmuldiagonal'", m,d);
    1044         350 :   for (j=1; j<lx; j++) gel(y,j) = RgC_Rg_mul(gel(m,j), gel(d,j));
    1045          70 :   return y;
    1046             : }
    1047             : 
    1048             : /* compute A*B assuming the result is a diagonal matrix */
    1049             : GEN
    1050           7 : matmultodiagonal(GEN A, GEN B)
    1051             : {
    1052           7 :   long i, j, hA, hB, lA = lg(A), lB = lg(B);
    1053           7 :   GEN y = matid(lB-1);
    1054             : 
    1055           7 :   if (typ(A) != t_MAT) pari_err_TYPE("matmultodiagonal",A);
    1056           7 :   if (typ(B) != t_MAT) pari_err_TYPE("matmultodiagonal",B);
    1057           7 :   hA = (lA == 1)? lB: lgcols(A);
    1058           7 :   hB = (lB == 1)? lA: lgcols(B);
    1059           7 :   if (lA != hB || lB != hA) pari_err_OP("operation 'matmultodiagonal'", A,B);
    1060          56 :   for (i=1; i<lB; i++)
    1061             :   {
    1062          49 :     GEN z = gen_0;
    1063         392 :     for (j=1; j<lA; j++) z = gadd(z, gmul(gcoeff(A,i,j),gcoeff(B,j,i)));
    1064          49 :     gcoeff(y,i,i) = z;
    1065             :   }
    1066           7 :   return y;
    1067             : }
    1068             : 
    1069             : /* [m[1,1], ..., m[l,l]], internal */
    1070             : GEN
    1071      885669 : RgM_diagonal_shallow(GEN m)
    1072             : {
    1073      885669 :   long i, lx = lg(m);
    1074      885669 :   GEN y = cgetg(lx,t_VEC);
    1075     2909245 :   for (i=1; i<lx; i++) gel(y, i) = gcoeff(m,i,i);
    1076      885708 :   return y;
    1077             : }
    1078             : 
    1079             : /* same, public function */
    1080             : GEN
    1081           0 : RgM_diagonal(GEN m)
    1082             : {
    1083           0 :   long i, lx = lg(m);
    1084           0 :   GEN y = cgetg(lx,t_VEC);
    1085           0 :   for (i=1; i<lx; i++) gel(y,i) = gcopy(gcoeff(m,i,i));
    1086           0 :   return y;
    1087             : }

Generated by: LCOV version 1.16