source: git/Singular/mpr_base.cc @ fdc537

spielwiese
Last change on this file since fdc537 was 27b040, checked in by Hans Schönemann <hannes@…>, 24 years ago
*** empty log message *** git-svn-id: file:///usr/local/Singular/svn/trunk@4170 2c84dea3-7e68-4137-9b89-c4e89433aadc
  • Property mode set to 100644
File size: 73.7 KB
Line 
1/****************************************
2*  Computer Algebra System SINGULAR     *
3****************************************/
4/* $Id: mpr_base.cc,v 1.20 2000-03-01 14:00:05 Singular Exp $ */
5
6/*
7 * ABSTRACT - multipolynomial resultants - resultant matrices
8 *            ( sparse, dense, u-resultant solver )
9 */
10
11#include <math.h>
12#include <limits.h>
13
14#include "mod2.h"
15
16//-> includes
17#include "structs.h"
18#include "polys.h"
19#include "ideals.h"
20#include "ipid.h"
21#include "ipshell.h"
22#include "febase.h"
23#include "mmemory.h"
24#include "intvec.h"
25#include "matpol.h"
26#include "numbers.h"
27#include "longalg.h"
28#ifdef HAVE_FACTORY
29#include "clapsing.h"
30#endif
31#include "sparsmat.h"
32
33#include "mpr_global.h"
34#include "mpr_base.h"
35#include "mpr_numeric.h"
36#if HAVE_ASO == 10
37#include "mpr_base.aso"
38#else
39#define AllocSizeOf(X) Alloc(sizeof(X))
40#define FreeSizeOf(X,Y) Free(X,sizeof(Y))
41#endif
42//<-
43
44extern void nPrint(number n);  // for debugging output
45
46//%s
47//-----------------------------------------------------------------------------
48//-------------- sparse resultant matrix --------------------------------------
49//-----------------------------------------------------------------------------
50
51//-> definitions
52
53//#define mprTEST
54//#define mprMINKSUM
55
56#define MAXPOINTS      10000
57#define MAXINITELEMS   256
58#define LIFT_COOR      50000   // siRand() % LIFT_COOR gives random lift value
59#define SCALEDOWN      100.0  // lift value scale down for linear program
60#define MINVDIST       0.0
61#define RVMULT         0.0001 // multiplicator for random shift vector
62#define MAXRVVAL       50000
63#define MAXVARS        100
64//<-
65
66//-> sparse resultant matrix
67
68/* set of points */
69class pointSet;
70
71
72
73/* sparse resultant matrix class */
74class resMatrixSparse : virtual public resMatrixBase
75{
76public:
77  resMatrixSparse( const ideal _gls, const int special = SNONE );
78  ~resMatrixSparse();
79
80  // public interface according to base class resMatrixBase
81  const ideal getMatrix();
82
83  /** Fills in resMat[][] with evpoint[] and gets determinant
84   * uRPos[i][1]: row of matrix
85   * uRPos[i][idelem+1]: col of u(0)
86   *  uRPos[i][2..idelem]: col of u(1) .. u(n)
87   *  i= 1 .. numSet0
88   */
89  const number getDetAt( const number* evpoint );
90
91  const poly getUDet( const number* evpoint );
92
93private:
94  resMatrixSparse( const resMatrixSparse & );
95
96  void randomVector( const int dim, mprfloat shift[] );
97
98  /** Row Content Function
99   * Finds the largest i such that F[i] is a point, F[i]= a[ij] in A[i] for some j.
100   * Returns -1 iff the point vert does not lie in a cell
101   */
102  int RC( pointSet **pQ, pointSet *E, int vert, mprfloat shift[] );
103
104  /* Remaps a result of LP to the according point set Qi.
105   * Returns false iff remaping was not possible, otherwise true.
106   */
107  bool remapXiToPoint( const int indx, pointSet **pQ, int *set, int *vtx );
108
109  /** create coeff matrix
110   * uRPos[i][1]: row of matrix
111   * uRPos[i][idelem+1]: col of u(0)
112   *  uRPos[i][2..idelem]: col of u(1) .. u(n)
113   *  i= 1 .. numSet0
114   * Returns the dimension of the matrix or -1 in case of an error
115   */
116  int createMatrix( pointSet *E );
117
118  pointSet * minkSumAll( pointSet **pQ, int numq, int dim );
119  pointSet * minkSumTwo( pointSet *Q1, pointSet *Q2, int dim );
120
121private:
122  ideal gls;
123
124  int n, idelem;     // number of variables, polynoms
125  int numSet0;       // number of elements in S0
126  int msize;         // size of matrix
127
128  intvec *uRPos;
129
130  ideal rmat;        // sparse matrix representation
131
132  simplex * LP;      // linear programming stuff
133};
134//<-
135
136//-> typedefs and structs
137poly monomAt( poly p, int i );
138
139typedef unsigned int Coord_t;
140
141struct setID
142{
143  int set;
144  int pnt;
145};
146
147struct onePoint
148{
149  Coord_t * point;             // point[0] is unused, maxial dimension is MAXVARS+1
150  setID rc;                    // filled in by Row Content Function
151  struct onePoint * rcPnt;     // filled in by Row Content Function
152};
153
154typedef struct onePoint * onePointP;
155
156/* sparse matrix entry */
157struct _entry
158{
159  number num;
160  int col;
161  struct _entry * next;
162};
163
164typedef struct _entry * entry;
165//<-
166
167//-> class pointSet
168class pointSet
169{
170private:
171  onePointP *points;     // set of onePoint's, index [1..num], supports of monoms
172  bool lifted;
173
174public:
175  int num;               // number of elements in points
176  int max;               // maximal entries in points, i.e. allocated mem
177  int dim;               // dimension, i.e. valid coord entries in point
178  int index;             // should hold unique identifier of point set
179
180  pointSet( const int _dim, const int _index= 0, const int count= MAXINITELEMS );
181   ~pointSet();
182
183  // pointSet.points[i] equals pointSet[i]
184  inline const onePointP operator[] ( const int index );
185
186  /** Adds a point to pointSet, copy vert[0,...,dim] ot point[num+1][0,...,dim].
187   * Returns false, iff additional memory was allocated ( i.e. num >= max )
188   * else returns true
189   */
190  bool addPoint( const onePointP vert );
191
192  /** Adds a point to pointSet, copy vert[0,...,dim] ot point[num+1][0,...,dim].
193   * Returns false, iff additional memory was allocated ( i.e. num >= max )
194   * else returns true
195   */
196  bool addPoint( const Exponent_t * vert );
197
198  /** Adds a point to pointSet, copy vert[0,...,dim] ot point[num+1][0,...,dim].
199   * Returns false, iff additional memory was allocated ( i.e. num >= max )
200   * else returns true
201   */
202  bool addPoint( const Coord_t * vert );
203
204  /* Removes the point at intex indx */
205  bool removePoint( const int indx );
206
207  /** Adds point to pointSet, iff pointSet \cap point = \emptyset.
208   * Returns true, iff added, else false.
209   */
210  bool mergeWithExp( const onePointP vert );
211
212  /** Adds point to pointSet, iff pointSet \cap point = \emptyset.
213   * Returns true, iff added, else false.
214   */
215  bool mergeWithExp( const Exponent_t * vert );
216
217  /* Adds support of poly p to pointSet, iff pointSet \cap point = \emptyset. */
218  void mergeWithPoly( const poly p );
219
220  /* Returns the row polynom multiplicator in vert[] */
221  void getRowMP( const int indx, Exponent_t * vert );
222
223  /* Returns index of supp(LT(p)) in pointSet. */
224  int getExpPos( const poly p );
225
226  /** sort lex
227   */
228  void sort();
229
230  /** Lifts the point set using sufficiently generic linear lifting
231   * homogeneous forms l[1]..l[dim] in Z. Every l[i] is of the form
232   * L1x1+...+Lnxn, for generic L1..Ln in Z.
233   *
234   * Lifting raises dimension by one!
235   */
236  void lift( int *l= NULL );     // !! increments dim by 1
237  void unlift() { dim--; lifted= false; }
238
239private:
240  pointSet( const pointSet & );
241
242  /** points[a] < points[b] ? */
243  inline bool smaller( int, int );
244
245  /** points[a] > points[b] ? */
246  inline bool larger( int, int );
247
248  /** Checks, if more mem is needed ( i.e. num >= max ),
249   * returns false, if more mem was allocated, else true
250   */
251  inline bool checkMem();
252};
253//<-
254
255//-> class convexHull
256/* Compute convex hull of given exponent set */
257class convexHull
258{
259public:
260  convexHull( simplex * _pLP ) : pLP(_pLP) {}
261  ~convexHull() {}
262
263  /** Computes the point sets of the convex hulls of the supports given
264   * by the polynoms in gls.
265   * Returns Q[].
266   */
267  pointSet ** newtonPolytopesP( const ideal gls );
268  ideal newtonPolytopesI( const ideal gls );
269
270private:
271  /** Returns true iff the support of poly pointPoly is inside the
272   * convex hull of all points given by the  support of poly p.
273   */
274  bool inHull(poly p, poly pointPoly, int m, int site);
275
276private:
277  pointSet **Q;
278  int n;
279  simplex * pLP;
280};
281//<-
282
283//-> class mayanPyramidAlg
284/* Compute all lattice points in a given convex hull */
285class mayanPyramidAlg
286{
287public:
288  mayanPyramidAlg( simplex * _pLP ) : n(pVariables), pLP(_pLP) {}
289  ~mayanPyramidAlg() {}
290
291  /** Drive Mayan Pyramid Algorithm.
292   * The Alg computes conv(Qi[]+shift[]).
293   */
294  pointSet * getInnerPoints( pointSet **_Qi, mprfloat _shift[] );
295
296private:
297
298  /** Recursive Mayan Pyramid algorithm for directly computing MinkowskiSum
299   * lattice points for (n+1)-fold MinkowskiSum of given point sets Qi[].
300   * Recursively for range of dim: dim in [0..n); acoords[0..var) fixed.
301   * Stores only MinkowskiSum points of udist > 0: done by storeMinkowskiSumPoints.
302   */
303  void runMayanPyramid( int dim );
304
305  /**  Compute v-distance via Linear Programing
306   * Linear Program finds the v-distance of the point in accords[].
307   * The v-distance is the distance along the direction v to boundary of
308   * Minkowski Sum of Qi (here vector v is represented by shift[]).
309   * Returns the v-distance or -1.0 if an error occured.
310   */
311  mprfloat vDistance( Coord_t * acoords, int dim );
312
313  /** LP for finding min/max coord in MinkowskiSum, given previous coors.
314   * Assume MinkowskiSum in non-negative quadrants
315   * coor in [0,n); fixed coords in acoords[0..coor)
316   */
317  void mn_mx_MinkowskiSum( int dim, Coord_t *minR, Coord_t *maxR );
318
319  /**  Stores point in E->points[pt], iff v-distance != 0
320   * Returns true iff point was stored, else flase
321   */
322  bool storeMinkowskiSumPoint();
323
324private:
325  pointSet **Qi;
326  pointSet *E;
327  mprfloat *shift;
328
329  int n,idelem;
330
331  Coord_t acoords[MAXVARS+2];
332
333  simplex * pLP;
334};
335//<-
336
337#ifndef ASO_GENERATE
338
339//-> debug output stuff
340#if defined(mprDEBUG_PROT) || defined(mprDEBUG_ALL)
341void print_mat(mprfloat **a, int maxrow, int maxcol)
342{
343  int i, j;
344
345  for (i = 1; i <= maxrow; i++)
346  {
347    Print("[");
348    for (j = 1; j <= maxcol; j++) Print("% 7.2f, ", a[i][j]);
349    Print("],\n");
350  }
351}
352void print_bmat(mprfloat **a, int nrows, int ncols, int N, int *iposv)
353{
354  int i, j;
355
356  printf("Output matrix from LinProg");
357  for (i = 1; i <= nrows; i++)
358  {
359    printf("\n[ ");
360    if (i == 1) printf("  ");
361    else if (iposv[i-1] <= N) printf("X%d", iposv[i-1]);
362    else printf("Y%d", iposv[i-1]-N+1);
363    for (j = 1; j <= ncols; j++) printf(" %7.2f ",(double)a[i][j]);
364    printf(" ]");
365  } printf("\n");
366  fflush(stdout);
367}
368
369void print_exp( const onePointP vert, int n )
370{
371  int i;
372  for ( i= 1; i <= n; i++ )
373  {
374    Print(" %d",vert->point[i] );
375#ifdef LONG_OUTPUT
376    if ( i < n ) Print(", ");
377#endif
378  }
379}
380void print_matrix( matrix omat )
381{
382  int i,j;
383  int val;
384  Print(" matrix m[%d][%d]=(\n",MATROWS( omat ),MATCOLS( omat ));
385  for ( i= 1; i <= MATROWS( omat ); i++ )
386  {
387    for ( j= 1; j <= MATCOLS( omat ); j++ )
388    {
389      if ( (MATELEM( omat, i, j)!=NULL)
390      && (!nIsZero(pGetCoeff( MATELEM( omat, i, j)))))
391      {
392        val= nInt(pGetCoeff( MATELEM( omat, i, j) ));
393        if ( i==MATROWS(omat) && j==MATCOLS(omat) )
394        {
395          Print("%d ",val);
396        }
397        else
398        {
399          Print("%d, ",val);
400        }
401      }
402      else
403      {
404        if ( i==MATROWS(omat) && j==MATCOLS(omat) )
405        {
406          PrintS("  0");
407        }
408        else
409        {
410          PrintS("  0, ");
411        }
412      }
413    }
414    PrintLn();
415  }
416  PrintS(");\n");
417}
418#endif
419//<-
420
421//-> pointSet::*
422pointSet::pointSet( const int _dim, const int _index, const int count )
423  : num(0), max(count), dim(_dim), index(_index)
424{
425  int i;
426  points = (onePointP *)Alloc( (count+1) * sizeof(onePointP) );
427  for ( i= 0; i <= max; i++ )
428  {
429    points[i]= (onePointP)Alloc( sizeof(onePoint) );
430    points[i]->point= (Coord_t *)Alloc0( (dim+2) * sizeof(Coord_t) );
431  }
432  lifted= false;
433}
434
435pointSet::~pointSet()
436{
437  int i;
438  int fdim= lifted ? dim+1 : dim+2;
439  for ( i= 0; i <= max; i++ )
440  {
441    Free( (ADDRESS) points[i]->point, fdim * sizeof(Coord_t) );
442    Free( (ADDRESS) points[i], sizeof(onePoint) );
443  }
444  Free( (ADDRESS) points, (max+1) * sizeof(onePointP) );
445}
446
447inline const onePointP pointSet::operator[] ( const int index )
448{
449  assume( index > 0 && index <= num );
450  return points[index];
451}
452
453inline bool pointSet::checkMem()
454{
455  if ( num >= max )
456  {
457    int i;
458    int fdim= lifted ? dim+1 : dim+2;
459    points= (onePointP*)ReAlloc( points,
460                                 (max+1) * sizeof(onePointP),
461                                 (2*max + 1) * sizeof(onePointP) );
462    for ( i= max+1; i <= max*2; i++ )
463    {
464      points[i]= (onePointP)Alloc( sizeof(struct onePoint) );
465      points[i]->point= (Coord_t *)Alloc0( fdim * sizeof(Coord_t) );
466    }
467    max*= 2;
468    mprSTICKYPROT(ST_SPARSE_MEM);
469    return false;
470  }
471  return true;
472}
473
474bool pointSet::addPoint( const onePointP vert )
475{
476  int i;
477  bool ret;
478  num++;
479  ret= checkMem();
480  points[num]->rcPnt= NULL;
481  for ( i= 1; i <= dim; i++ ) points[num]->point[i]= vert->point[i];
482  return ret;
483}
484
485bool pointSet::addPoint( const Exponent_t * vert )
486{
487  int i;
488  bool ret;
489  num++;
490  ret= checkMem();
491  points[num]->rcPnt= NULL;
492  for ( i= 1; i <= dim; i++ ) points[num]->point[i]= (Coord_t) vert[i];
493  return ret;
494}
495
496bool pointSet::addPoint( const Coord_t * vert )
497{
498  int i;
499  bool ret;
500  num++;
501  ret= checkMem();
502  points[num]->rcPnt= NULL;
503  for ( i= 0; i < dim; i++ ) points[num]->point[i+1]= vert[i];
504  return ret;
505}
506
507bool pointSet::removePoint( const int indx )
508{
509  assume( indx > 0 && indx <= num );
510  if ( indx != num )
511  {
512    onePointP tmp;
513    tmp= points[indx];
514    points[indx]= points[num];
515    points[num]= tmp;
516  }
517  num--;
518
519  return true;
520}
521
522bool pointSet::mergeWithExp( const onePointP vert )
523{
524  int i,j;
525
526  for ( i= 1; i <= num; i++ )
527  {
528    for ( j= 1; j <= dim; j++ )
529      if ( points[i]->point[j] != vert->point[j] ) break;
530    if ( j > dim ) break;
531  }
532
533  if ( i > num )
534  {
535    addPoint( vert );
536    return true;
537  }
538  return false;
539}
540
541bool pointSet::mergeWithExp( const Exponent_t * vert )
542{
543  int i,j;
544
545  for ( i= 1; i <= num; i++ )
546  {
547    for ( j= 1; j <= dim; j++ )
548      if ( points[i]->point[j] != (Coord_t) vert[j] ) break;
549    if ( j > dim ) break;
550  }
551
552  if ( i > num )
553  {
554    addPoint( vert );
555    return true;
556  }
557  return false;
558}
559
560void pointSet::mergeWithPoly( const poly p )
561{
562  int i,j;
563  poly piter= p;
564  Exponent_t * vert;
565  vert= (Exponent_t *)Alloc( (dim+1) * sizeof(Exponent_t) );
566
567  while ( piter )
568  {
569    pGetExpV( piter, vert );
570
571    for ( i= 1; i <= num; i++ )
572    {
573      for ( j= 1; j <= dim; j++ )
574        if ( points[i]->point[j] != (Coord_t) vert[j] ) break;
575      if ( j > dim ) break;
576    }
577
578    if ( i > num )
579    {
580      addPoint( vert );
581    }
582
583    pIter( piter );
584  }
585  Free( (ADDRESS) vert, (dim+1) * sizeof(Exponent_t) );
586}
587
588int pointSet::getExpPos( const poly p )
589{
590  Exponent_t * vert;
591  int i,j;
592
593  // hier unschoen...
594  vert= (Exponent_t *)Alloc( (dim+1) * sizeof(Exponent_t) );
595
596  pGetExpV( p, vert );
597  for ( i= 1; i <= num; i++ )
598  {
599    for ( j= 1; j <= dim; j++ )
600      if ( points[i]->point[j] != (Coord_t) vert[j] ) break;
601    if ( j > dim ) break;
602  }
603  Free( (ADDRESS) vert, (dim+1) * sizeof(Exponent_t) );
604
605  if ( i > num ) return 0;
606  else return i;
607}
608
609void pointSet::getRowMP( const int indx, Exponent_t * vert )
610{
611  assume( indx > 0 && indx <= num && points[indx]->rcPnt );
612  int i;
613
614  vert[0]= 0;
615  for ( i= 1; i <= dim; i++ )
616    vert[i]= (Exponent_t)(points[indx]->point[i] - points[indx]->rcPnt->point[i]);
617}
618
619inline bool pointSet::smaller( int a, int b )
620{
621  int i;
622
623  for ( i= 1; i <= dim; i++ )
624  {
625    if ( points[a]->point[i] > points[b]->point[i] )
626    {
627      return false;
628    }
629    if ( points[a]->point[i] < points[b]->point[i] )
630    {
631      return true;
632    }
633  }
634
635 return false; // they are equal
636}
637
638inline bool pointSet::larger( int a, int b )
639{
640  int i;
641
642  for ( i= 1; i <= dim; i++ )
643  {
644    if ( points[a]->point[i] < points[b]->point[i] )
645    {
646      return false;
647    }
648    if ( points[a]->point[i] > points[b]->point[i] )
649    {
650      return true;
651    }
652  }
653
654 return false; // they are equal
655}
656
657void pointSet::sort()
658{
659  int i,j;
660  bool found= true;
661  onePointP tmp;
662
663  while ( found )
664  {
665    found= false;
666    for ( i= 1; i < num; i++ )
667    {
668      if ( larger( i, i+1 ) )
669      {
670        tmp= points[i];
671        points[i]= points[i+1];
672        points[i+1]= tmp;
673
674        found= true;
675      }
676    }
677  }
678}
679
680void pointSet::lift( int l[] )
681{
682  bool outerL= true;
683  int i, j;
684  int sum;
685
686  dim++;
687
688  if ( l==NULL )
689  {
690    outerL= false;
691    l= (int *)Alloc( (dim+1) * sizeof(int) ); // [1..dim-1]
692
693    for(i = 1; i < dim; i++)
694    {
695      l[i]= 1 + siRand() % LIFT_COOR;
696    }
697  }
698  for ( j=1; j <= num; j++ )
699  {
700    sum= 0;
701    for ( i=1; i < dim; i++ )
702    {
703      sum += (int)points[j]->point[i] * l[i];
704    }
705    points[j]->point[dim]= sum;
706  }
707
708#ifdef mprDEBUG_ALL
709  Print(" lift vector: ");
710  for ( j=1; j < dim; j++ ) Print(" %d ",l[j] );
711  PrintLn();
712#ifdef mprDEBUG_ALL
713  Print(" lifted points: \n");
714  for ( j=1; j <= num; j++ )
715  {
716    Print("%d: <",j);print_exp(points[j],dim);Print(">\n");
717  }
718  PrintLn();
719#endif
720#endif
721
722  lifted= true;
723
724  if ( !outerL ) Free( (ADDRESS) l, (dim+1) * sizeof(int) );
725}
726//<-
727
728//-> global functions
729// Returns the monom at pos i in poly p
730poly monomAt( poly p, int i )
731{
732  assume( i > 0 );
733  poly iter= p;
734  for ( int j= 1; (j < i) && iter; j++ ) iter= pIter(iter);
735  return iter;
736}
737//<-
738
739//-> convexHull::*
740bool convexHull::inHull(poly p, poly pointPoly, int m, int site)
741{
742  int i, j, col;
743
744  pLP->m = n+1;
745  pLP->n = m;                // this includes col of cts
746
747  pLP->LiPM[1][1] = +0.0;
748  pLP->LiPM[1][2] = +1.0;        // optimize (arbitrary) var
749  pLP->LiPM[2][1] = +1.0;
750  pLP->LiPM[2][2] = -1.0;         // lambda vars sum up to 1
751
752  for ( j=3; j <= pLP->n; j++)
753  {
754    pLP->LiPM[1][j] = +0.0;
755    pLP->LiPM[2][j] = -1.0;
756  }
757
758  for( i= 1; i <= n; i++) {        // each row constraints one coor
759    pLP->LiPM[i+2][1] = (mprfloat)pGetExp(pointPoly,i);
760    col = 2;
761    for( j= 1; j <= m; j++ )
762    {
763      if( j != site )
764      {
765        pLP->LiPM[i+2][col] = -(mprfloat)pGetExp( monomAt(p,j), i );
766        col++;
767      }
768    }
769  }
770
771#ifdef mprDEBUG_ALL
772  PrintS("Matrix of Linear Programming\n");
773  print_mat( pLP->LiPM, pLP->m+1,pLP->n);
774#endif
775
776  pLP->m3= pLP->m;
777
778  pLP->compute();
779
780  return (pLP->icase == 0);
781}
782
783// mprSTICKYPROT:
784// ST_SPARSE_VADD: new vertex of convex hull added
785// ST_SPARSE_VREJ: point rejected (-> inside hull)
786pointSet ** convexHull::newtonPolytopesP( const ideal gls )
787{
788  int i, j, k;
789  int m;  // Anzahl der Exponentvektoren im i-ten Polynom (gls->m)[i] des Ideals gls
790  int idelem= IDELEMS(gls);
791  Exponent_t * vert;
792
793  n= pVariables;
794  vert= (Exponent_t *)Alloc( (idelem+1) * sizeof(Exponent_t) );
795
796  Q = (pointSet **)Alloc( idelem * sizeof(pointSet*) );        // support hulls
797  for ( i= 0; i < idelem; i++ )
798    Q[i] = new pointSet( pVariables, i+1, pLength((gls->m)[i])+1 );
799
800  for( i= 0; i < idelem; i++ )
801  {
802    k=1;
803    m = pLength( (gls->m)[i] );
804
805    poly p= (gls->m)[i];
806    for( j= 1; j <= m; j++) {  // für jeden Exponentvektor
807      if( !inHull( (gls->m)[i], p, m, j ) )
808      {
809        pGetExpV( p, vert );
810        Q[i]->addPoint( vert );
811        k++;
812        mprSTICKYPROT(ST_SPARSE_VADD);
813      }
814      else
815      {
816        mprSTICKYPROT(ST_SPARSE_VREJ);
817      }
818      pIter( p );
819    } // j
820    mprSTICKYPROT("\n");
821  } // i
822
823  Free( (ADDRESS) vert, (idelem+1) * sizeof(Exponent_t) );
824
825#ifdef mprDEBUG_PROT
826  PrintLn();
827  for( i= 0; i < idelem; i++ )
828  {
829    Print(" \\Conv(Qi[%d]): #%d\n", i,Q[i]->num );
830    for ( j=1; j <= Q[i]->num; j++ )
831    {
832      Print("%d: <",j);print_exp( (*Q[i])[j] , pVariables );PrintS(">\n");
833    }
834    PrintLn();
835  }
836#endif
837
838  return Q;
839}
840
841// mprSTICKYPROT:
842// ST_SPARSE_VADD: new vertex of convex hull added
843// ST_SPARSE_VREJ: point rejected (-> inside hull)
844ideal convexHull::newtonPolytopesI( const ideal gls )
845{
846  int i, j;
847  int m;  // Anzahl der Exponentvektoren im i-ten Polynom (gls->m)[i] des Ideals gls
848  int idelem= IDELEMS(gls);
849  ideal id;
850  poly p,pid,pd;
851  Exponent_t * vert;
852
853  n= pVariables;
854  vert= (Exponent_t *)Alloc( (idelem+1) * sizeof(Exponent_t) );
855  id= idInit( idelem, 1 );
856
857  for( i= 0; i < idelem; i++ )
858  {
859    m = pLength( (gls->m)[i] );
860
861    p= (gls->m)[i];
862    for( j= 1; j <= m; j++) {  // für jeden Exponentvektor
863      if( !inHull( (gls->m)[i], p, m, j ) )
864      {
865        if ( (id->m)[i] == NULL )
866        {
867          (id->m)[i]= pHead(p);
868          pid=(id->m)[i];
869        }
870        else
871        {
872          pNext(pid)= pHead(p);
873          pIter(pid);
874          pNext(pid)= NULL;
875        }
876        mprSTICKYPROT(ST_SPARSE_VADD);
877      }
878      else
879      {
880        mprSTICKYPROT(ST_SPARSE_VREJ);
881      }
882      pIter( p );
883    } // j
884    mprSTICKYPROT("\n");
885  } // i
886
887  Free( (ADDRESS) vert, (idelem+1) * sizeof(Exponent_t) );
888
889#ifdef mprDEBUG_PROT
890  PrintLn();
891  for( i= 0; i < idelem; i++ )
892  {
893  }
894#endif
895
896  return id;
897}
898//<-
899
900//-> mayanPyramidAlg::*
901pointSet * mayanPyramidAlg::getInnerPoints( pointSet **_Qi, mprfloat _shift[] )
902{
903  int i;
904
905  Qi= _Qi;
906  shift= _shift;
907
908  E= new pointSet( Qi[0]->dim ); // E has same dim as Qi[...]
909
910  for ( i= 0; i < MAXVARS+2; i++ ) acoords[i]= 0;
911
912  runMayanPyramid(0);
913
914  mprSTICKYPROT("\n");
915
916  return E;
917}
918
919mprfloat mayanPyramidAlg::vDistance( Coord_t * acoords, int dim )
920{
921  int i, ii, j, k, col, r;
922  int numverts, cols;
923
924  numverts = 0;
925  for( i=0; i<=n; i++)
926  {
927    numverts += Qi[i]->num;
928  }
929  cols = numverts + 2;
930
931  //if( dim < 1 || dim > n )
932  //  WerrorS("mayanPyramidAlg::vDistance: Known coords dim off range");
933
934  pLP->LiPM[1][1] = 0.0;
935  pLP->LiPM[1][2] = 1.0;        // maximize
936  for( j=3; j<=cols; j++) pLP->LiPM[1][j] = 0.0;
937
938  for( i=0; i <= n; i++ )
939  {
940    pLP->LiPM[i+2][1] = 1.0;
941    pLP->LiPM[i+2][2] = 0.0;
942  }
943  for( i=1; i<=dim; i++)
944  {
945    pLP->LiPM[n+2+i][1] = (mprfloat)(acoords[i-1]);
946    pLP->LiPM[n+2+i][2] = -shift[i];
947  }
948
949  ii = -1;
950  col = 2;
951  for ( i= 0; i <= n; i++ )
952  {
953    ii++;
954    for( k= 1; k <= Qi[ii]->num; k++ )
955    {
956      col++;
957      for ( r= 0; r <= n; r++ )
958      {
959        if ( r == i ) pLP->LiPM[r+2][col] = -1.0;
960        else pLP->LiPM[r+2][col] = 0.0;
961      }
962      for( r= 1; r <= dim; r++ )
963        pLP->LiPM[r+n+2][col] = -(mprfloat)((*Qi[ii])[k]->point[r]);
964    }
965  }
966
967  if( col != cols)
968    Werror("mayanPyramidAlg::vDistance:"
969           "setting up matrix for udist: col %d != cols %d",col,cols);
970
971  pLP->m = n+dim+1;
972  pLP->m3= pLP->m;
973  pLP->n=cols-1;
974
975#ifdef mprDEBUG_ALL
976  Print("vDistance LP, known koords dim=%d, constr %d, cols %d, acoords= ",
977        dim,pLP->m,cols);
978  for( i= 0; i < dim; i++ )
979    Print(" %d",acoords[i]);
980  PrintLn();
981  print_mat( pLP->LiPM, pLP->m+1, cols);
982#endif
983
984  pLP->compute();
985
986#ifdef mprDEBUG_ALL
987  Print("LP returns matrix\n");
988  print_bmat( pLP->LiPM, pLP->m+1, cols+1-pLP->m, cols, pLP->iposv);
989#endif
990
991  if( pLP->icase != 0 )
992  {  // check for errors
993    WerrorS("mayanPyramidAlg::vDistance:");
994    if( pLP->icase == 1 )
995      WerrorS(" Unbounded v-distance: probably 1st v-coor=0");
996    else if( pLP->icase == -1 )
997      WerrorS(" Infeasible v-distance");
998    else
999      WerrorS(" Unknown error");
1000    return -1.0;
1001  }
1002
1003  return pLP->LiPM[1][1];
1004}
1005
1006void  mayanPyramidAlg::mn_mx_MinkowskiSum( int dim, Coord_t *minR, Coord_t *maxR )
1007{
1008  int i, j, k, cols, cons;
1009  int la_cons_row;
1010
1011  cons = n+dim+2;
1012
1013  // first, compute minimum
1014  //
1015
1016  // common part of the matrix
1017  pLP->LiPM[1][1] = 0.0;
1018  for( i=2; i<=n+2; i++)
1019  {
1020    pLP->LiPM[i][1] = 1.0;        // 1st col
1021    pLP->LiPM[i][2] = 0.0;        // 2nd col
1022  }
1023
1024  la_cons_row = 1;
1025  cols = 2;
1026  for( i=0; i<=n; i++)
1027  {
1028    la_cons_row++;
1029    for( j=1; j<= Qi[i]->num; j++)
1030    {
1031      cols++;
1032      pLP->LiPM[1][cols] = 0.0;        // set 1st row 0
1033      for( k=2; k<=n+2; k++)
1034      {  // lambdas sum up to 1
1035        if( k != la_cons_row) pLP->LiPM[k][cols] = 0.0;
1036        else pLP->LiPM[k][cols] = -1.0;
1037      }
1038      for( k=1; k<=n; k++)
1039        pLP->LiPM[k+n+2][cols] = -(mprfloat)((*Qi[i])[j]->point[k]);
1040    } // j
1041  } // i
1042
1043  for( i= 0; i < dim; i++ )
1044  {                // fixed coords
1045    pLP->LiPM[i+n+3][1] = acoords[i];
1046    pLP->LiPM[i+n+3][2] = 0.0;
1047  }
1048  pLP->LiPM[dim+n+3][1] = 0.0;
1049
1050
1051  pLP->LiPM[1][2] = -1.0;                        // minimize
1052  pLP->LiPM[dim+n+3][2] = 1.0;
1053
1054#ifdef mprDEBUG_ALL
1055  Print("\nThats the matrix for minR, dim= %d, acoords= ",dim);
1056  for( i= 0; i < dim; i++ )
1057    Print(" %d",acoords[i]);
1058  PrintLn();
1059  print_mat( pLP->LiPM, cons+1, cols);
1060#endif
1061
1062  // simplx finds MIN for obj.fnc, puts it in [1,1]
1063  pLP->m= cons;
1064  pLP->n= cols-1;
1065  pLP->m3= cons;
1066
1067  pLP->compute();
1068
1069  if ( pLP->icase != 0 )
1070  { // check for errors
1071    if( pLP->icase < 0)
1072      WerrorS(" mn_mx_MinkowskiSum: LinearProgram: minR: infeasible");
1073    else if( pLP->icase > 0)
1074      WerrorS(" mn_mx_MinkowskiSum: LinearProgram: minR: unbounded");
1075  }
1076
1077  *minR = (Coord_t)( -pLP->LiPM[1][1] + 1.0 - SIMPLEX_EPS );
1078
1079  // now compute maximum
1080  //
1081
1082  // common part of the matrix again
1083  pLP->LiPM[1][1] = 0.0;
1084  for( i=2; i<=n+2; i++)
1085  {
1086    pLP->LiPM[i][1] = 1.0;
1087    pLP->LiPM[i][2] = 0.0;
1088  }
1089  la_cons_row = 1;
1090  cols = 2;
1091  for( i=0; i<=n; i++)
1092  {
1093    la_cons_row++;
1094    for( j=1; j<=Qi[i]->num; j++)
1095    {
1096      cols++;
1097      pLP->LiPM[1][cols] = 0.0;
1098      for( k=2; k<=n+2; k++)
1099      {
1100        if( k != la_cons_row) pLP->LiPM[k][cols] = 0.0;
1101        else pLP->LiPM[k][cols] = -1.0;
1102      }
1103      for( k=1; k<=n; k++)
1104        pLP->LiPM[k+n+2][cols] = -(mprfloat)((*Qi[i])[j]->point[k]);
1105    } // j
1106  }  // i
1107
1108  for( i= 0; i < dim; i++ )
1109  {                // fixed coords
1110    pLP->LiPM[i+n+3][1] = acoords[i];
1111    pLP->LiPM[i+n+3][2] = 0.0;
1112  }
1113  pLP->LiPM[dim+n+3][1] = 0.0;
1114
1115  pLP->LiPM[1][2] = 1.0;                      // maximize
1116  pLP->LiPM[dim+n+3][2] = 1.0;                // var = sum of pnt coords
1117
1118#ifdef mprDEBUG_ALL
1119  Print("\nThats the matrix for maxR, dim= %d\n",dim);
1120  print_mat( pLP->LiPM, cons+1, cols);
1121#endif
1122
1123  pLP->m= cons;
1124  pLP->n= cols-1;
1125  pLP->m3= cons;
1126
1127  // simplx finds MAX for obj.fnc, puts it in [1,1]
1128  pLP->compute();
1129
1130  if ( pLP->icase != 0 )
1131  {
1132    if( pLP->icase < 0)
1133      WerrorS(" mn_mx_MinkowskiSum: LinearProgram: maxR: infeasible");
1134    else if( pLP->icase > 0)
1135      WerrorS(" mn_mx_MinkowskiSum: LinearProgram: maxR: unbounded");
1136  }
1137
1138  *maxR = (Coord_t)( pLP->LiPM[1][1] + SIMPLEX_EPS );
1139
1140#ifdef mprDEBUG_ALL
1141  Print("  Range for dim=%d: [%d,%d]\n", dim, *minR, *maxR);
1142#endif
1143}
1144
1145// mprSTICKYPROT:
1146// ST_SPARSE_VREJ: rejected point
1147// ST_SPARSE_VADD: added point to set
1148bool mayanPyramidAlg::storeMinkowskiSumPoint()
1149{
1150  mprfloat dist;
1151
1152  // determine v-distance of point pt
1153  dist= vDistance( &(acoords[0]), n );
1154
1155  // store only points with v-distance > minVdist
1156  if( dist <= MINVDIST + SIMPLEX_EPS )
1157  {
1158    mprSTICKYPROT(ST_SPARSE_VREJ);
1159    return false;
1160  }
1161
1162  E->addPoint( &(acoords[0]) );
1163  mprSTICKYPROT(ST_SPARSE_VADD);
1164
1165  return true;
1166}
1167
1168// mprSTICKYPROT:
1169// ST_SPARSE_MREC1: recurse
1170// ST_SPARSE_MREC2: recurse with extra points
1171// ST_SPARSE_MPEND: end
1172void mayanPyramidAlg::runMayanPyramid( int dim )
1173{
1174  Coord_t minR, maxR;
1175  mprfloat dist;
1176
1177  // step 3
1178  mn_mx_MinkowskiSum( dim, &minR, &maxR );
1179
1180#ifdef mprDEBUG_ALL
1181  int i;
1182  for( i=0; i <= dim; i++) Print("acoords[%d]=%d ",i,(int)acoords[i]);
1183  Print(":: [%d,%d]\n", minR, maxR);
1184#endif
1185
1186  // step 5 -> terminate
1187  if( dim == n-1 )
1188  {
1189    int lastKilled = 0;
1190    // insert points
1191    acoords[dim] = minR;
1192    while( acoords[dim] <= maxR )
1193    {
1194      if( !storeMinkowskiSumPoint() )
1195        lastKilled++;
1196      acoords[dim]++;
1197    }
1198    mprSTICKYPROT(ST_SPARSE_MPEND);
1199    return;
1200  }
1201
1202  // step 4 -> recurse at step 3
1203  acoords[dim] = minR;
1204  while ( acoords[dim] <= maxR )
1205  {
1206    if ( (acoords[dim] > minR) && (acoords[dim] <= maxR) )
1207    {     // acoords[dim] >= minR  ??
1208      mprSTICKYPROT(ST_SPARSE_MREC1);
1209      runMayanPyramid( dim + 1 );         // recurse with higer dimension
1210    }
1211    else
1212    {
1213      // get v-distance of pt
1214      dist= vDistance( &(acoords[0]), dim + 1 );// dim+1 == known coordinates
1215
1216      if( dist >= SIMPLEX_EPS )
1217      {
1218        mprSTICKYPROT(ST_SPARSE_MREC2);
1219        runMayanPyramid( dim + 1 );       // recurse with higer dimension
1220      }
1221    }
1222    acoords[dim]++;
1223  } // while
1224}
1225//<-
1226
1227//-> resMatrixSparse::*
1228bool resMatrixSparse::remapXiToPoint( const int indx, pointSet **pQ, int *set, int *pnt )
1229{
1230  int i,n= pVariables;
1231  int loffset= 0;
1232  for ( i= 0; i <= n; i++ )
1233  {
1234    if ( (loffset < indx) && (indx <= pQ[i]->num + loffset) )
1235    {
1236      *set= i;
1237      *pnt= indx-loffset;
1238      return true;
1239    }
1240    else loffset+= pQ[i]->num;
1241  }
1242  return false;
1243}
1244
1245// mprSTICKYPROT
1246// ST_SPARSE_RC: point added
1247int resMatrixSparse::RC( pointSet **pQ, pointSet *E, int vert, mprfloat shift[] )
1248{
1249  int i, j, k,c ;
1250  int size;
1251  bool found= true;
1252  mprfloat cd;
1253  int set,pnt,onum;
1254  int bucket[MAXVARS+2];
1255  setID *optSum;
1256
1257  LP->n = 1;
1258  LP->m = n + n + 1;   // number of constrains
1259
1260  // fill in LP matrix
1261  for ( i= 0; i <= n; i++ )
1262  {
1263    size= pQ[i]->num;
1264    for ( k= 1; k <= size; k++ )
1265    {
1266      LP->n++;
1267
1268      // objective funtion, minimize
1269      LP->LiPM[1][LP->n] = - ( (mprfloat) (*pQ[i])[k]->point[pQ[i]->dim] / SCALEDOWN );
1270
1271      // lambdas sum up to 1
1272      for ( j = 0; j <= n; j++ )
1273      {
1274        if ( i==j )
1275          LP->LiPM[j+2][LP->n] = -1.0;
1276        else
1277          LP->LiPM[j+2][LP->n] = 0.0;
1278      }
1279
1280      // the points
1281      for ( j = 1; j <= n; j++ )
1282      {
1283        LP->LiPM[j+n+2][LP->n] =  - ( (mprfloat) (*pQ[i])[k]->point[j] );
1284      }
1285    }
1286  }
1287
1288  for ( j = 0; j <= n; j++ ) LP->LiPM[j+2][1] = 1.0;
1289  for ( j= 1; j <= n; j++ )
1290  {
1291    LP->LiPM[j+n+2][1]= (mprfloat)(*E)[vert]->point[j] - shift[j];
1292  }
1293  LP->n--;
1294
1295  LP->LiPM[1][1] = 0.0;
1296
1297#ifdef mprDEBUG_ALL
1298  PrintLn();
1299  Print(" n= %d, LP->m=M= %d, LP->n=N= %d\n",n,LP->m,LP->n);
1300  print_mat(LP->LiPM, LP->m+1, LP->n+1);
1301#endif
1302
1303  LP->m3= LP->m;
1304
1305  LP->compute();
1306
1307  if ( LP->icase < 0 )
1308  {
1309    // infeasibility: the point does not lie in a cell -> remove it
1310    return -1;
1311  }
1312
1313  // store result
1314  (*E)[vert]->point[E->dim]= (int)(-LP->LiPM[1][1] * SCALEDOWN);
1315
1316#ifdef mprDEBUG_ALL
1317  Print(" simplx returned %d, Objective value = %f\n", LP->icase, LP->LiPM[1][1]);
1318  //print_bmat(LP->LiPM, NumCons + 1, LP->n+1-NumCons, LP->n+1, LP->iposv); // ( rows= M+1, cols= N+1-m3 )
1319  //print_mat(LP->LiPM, NumCons+1, LP->n);
1320#endif
1321
1322#if 1
1323  // sort LP results
1324  while (found)
1325  {
1326    found=false;
1327    for ( i= 1; i < LP->m; i++ )
1328    {
1329      if ( LP->iposv[i] > LP->iposv[i+1] )
1330      {
1331
1332        c= LP->iposv[i];
1333        LP->iposv[i]=LP->iposv[i+1];
1334        LP->iposv[i+1]=c;
1335
1336        cd=LP->LiPM[i+1][1];
1337        LP->LiPM[i+1][1]=LP->LiPM[i+2][1];
1338        LP->LiPM[i+2][1]=cd;
1339
1340        found= true;
1341      }
1342    }
1343  }
1344#endif
1345
1346#ifdef mprDEBUG_ALL
1347  print_bmat(LP->LiPM, LP->m + 1, LP->n+1-LP->m, LP->n+1, LP->iposv);
1348  PrintS(" now split into sets\n");
1349#endif
1350
1351
1352  // init bucket
1353  for ( i= 0; i <= E->dim; i++ ) bucket[i]= 0;
1354  // remap results of LP to sets Qi
1355  c=0;
1356  optSum= (setID*)Alloc( (LP->m) * sizeof(struct setID) );
1357  for ( i= 0; i < LP->m; i++ )
1358  {
1359    //Print("% .15f\n",LP->LiPM[i+2][1]);
1360    if ( LP->LiPM[i+2][1] > 1e-12 )
1361    {
1362      if ( !remapXiToPoint( LP->iposv[i+1], pQ, &(optSum[c].set), &(optSum[c].pnt) ) )
1363      {
1364        Werror(" resMatrixSparse::RC: Found bad solution in LP: %d!",LP->iposv[i+1]);
1365        WerrorS(" resMatrixSparse::RC: remapXiToPoint faild!");
1366        return -1;
1367      }
1368      bucket[optSum[c].set]++;
1369      c++;
1370    }
1371  }
1372
1373  onum= c;
1374  // find last min in bucket[]: maximum i such that Fi is a point
1375  c= 0;
1376  for ( i= 1; i < E->dim; i++ )
1377  {
1378    if ( bucket[c] >= bucket[i] )
1379    {
1380      c= i;
1381    }
1382  }
1383  // find matching point set
1384  for ( i= onum - 1; i >= 0; i-- )
1385  {
1386    if ( optSum[i].set == c )
1387      break;
1388  }
1389  // store
1390  (*E)[vert]->rc.set= c;
1391  (*E)[vert]->rc.pnt= optSum[i].pnt;
1392  (*E)[vert]->rcPnt= (*pQ[c])[optSum[i].pnt];
1393  // count
1394  if ( (*E)[vert]->rc.set == linPolyS ) numSet0++;
1395
1396#ifdef mprDEBUG_PROT
1397  Print("\n Point E[%d] was <",vert);print_exp((*E)[vert],E->dim-1);Print(">, bucket={");
1398  for ( j= 0; j < E->dim; j++ )
1399  {
1400    Print(" %d",bucket[j]);
1401  }
1402  PrintS(" }\n optimal Sum: Qi ");
1403  for ( j= 0; j < LP->m; j++ )
1404  {
1405    Print(" [ %d, %d ]",optSum[j].set,optSum[j].pnt);
1406  }
1407  Print(" -> i= %d, j = %d\n",(*E)[vert]->rc.set,optSum[i].pnt);
1408#endif
1409
1410  // clean up
1411  Free( (ADDRESS) optSum, (LP->m) * sizeof(struct setID) );
1412
1413  mprSTICKYPROT(ST_SPARSE_RC);
1414
1415  return (int)(-LP->LiPM[1][1] * SCALEDOWN);
1416}
1417
1418// create coeff matrix
1419int resMatrixSparse::createMatrix( pointSet *E )
1420{
1421  // sparse matrix
1422  //    uRPos[i][1]: row of matrix
1423  //    uRPos[i][idelem+1]: col of u(0)
1424  //    uRPos[i][2..idelem]: col of u(1) .. u(n)
1425  //    i= 1 .. numSet0
1426  int i,j,epos;
1427  int rp,cp;
1428  poly rowp,epp;
1429  poly iterp;
1430  Exponent_t *epp_mon, *eexp;
1431
1432  epp_mon= (Exponent_t *)Alloc( (n+2) * sizeof(Exponent_t) );
1433  eexp= (Exponent_t *)Alloc0((pVariables+1)*sizeof(Exponent_t));
1434
1435  totDeg= numSet0;
1436
1437  mprSTICKYPROT2(" size of matrix: %d\n", E->num);
1438  mprSTICKYPROT2("  resultant deg: %d\n", numSet0);
1439
1440  uRPos= NewIntvec3( numSet0, pLength((gls->m)[0])+1, 0 );
1441
1442  // sparse Matrix represented as a module where
1443  // each poly is column vector ( pSetComp(p,k) gives the row )
1444  rmat= idInit( E->num, E->num );    // cols, rank= number of rows
1445  msize= E->num;
1446
1447  rp= 1;
1448  rowp= NULL;
1449  epp= pOne();
1450  for ( i= 1; i <= E->num; i++ )
1451  {       // for every row
1452    E->getRowMP( i, epp_mon );           // compute (p-a[ij]), (i,j) = RC(p)
1453    pSetExpV( epp, epp_mon );
1454
1455    //
1456    rowp= pMult( pCopy(epp), pCopy((gls->m)[(*E)[i]->rc.set]) );  // x^(p-a[ij]) * f(i)
1457
1458    cp= 2;
1459    // get column for every monomial in rowp and store it
1460    iterp= rowp;
1461    while ( iterp!=NULL )
1462    {
1463      epos= E->getExpPos( iterp );
1464      if ( epos == 0 )
1465      {
1466        // this can happen, if the shift vektor or the lift funktions
1467        // are not generically choosen.
1468        Werror("resMatrixSparse::createMatrix: Found exponent not in E, id %d, set [%d, %d]!",
1469               i,(*E)[i]->rc.set,(*E)[i]->rc.pnt);
1470        return i;
1471      }
1472      pSetExpV(iterp,eexp);
1473      pSetComp(iterp, epos );
1474      pSetm(iterp);
1475      if ( (*E)[i]->rc.set == linPolyS )
1476      { // store coeff positions
1477        IMATELEM(*uRPos,rp,cp)= epos;
1478        cp++;
1479      }
1480      pIter( iterp );
1481    } // while
1482    if ( (*E)[i]->rc.set == linPolyS )
1483    {   // store row
1484      IMATELEM(*uRPos,rp,1)= i-1;
1485      rp++;
1486    }
1487    (rmat->m)[i-1]= rowp;
1488  } // for
1489
1490  pDelete( &epp );
1491  Free( (ADDRESS) epp_mon, (n+2) * sizeof(Exponent_t) );
1492  Free( (ADDRESS) eexp, (pVariables+1)*sizeof(Exponent_t));
1493
1494#ifdef mprDEBUG_ALL
1495  if ( E->num <= 40 )
1496  {
1497    matrix mout= idModule2Matrix( idCopy(rmat) );
1498    print_matrix(mout);
1499  }
1500  for ( i= 1; i <= numSet0; i++ )
1501  {
1502    Print(" row  %d contains coeffs of f_%d\n",IMATELEM(*uRPos,i,1),linPolyS);
1503  }
1504  PrintS(" Sparse Matrix done\n");
1505#endif
1506
1507  return E->num;
1508}
1509
1510// find a sufficiently generic and small vector
1511void resMatrixSparse::randomVector( const int dim, mprfloat shift[] )
1512{
1513  int i,j;
1514  i= 1;
1515  time_t *tp = NULL;
1516
1517  while ( i <= dim )
1518  {
1519    shift[i]= (mprfloat) (RVMULT*(siRand()%MAXRVVAL)/(mprfloat)MAXRVVAL);
1520    i++;
1521    for ( j= 1; j < i-1; j++ )
1522    {
1523      if ( (shift[j] < shift[i-1] + SIMPLEX_EPS) && (shift[j] > shift[i-1] - SIMPLEX_EPS) )
1524      {
1525        i--;
1526        break;
1527      }
1528    }
1529  }
1530}
1531
1532pointSet * resMatrixSparse::minkSumTwo( pointSet *Q1, pointSet *Q2, int dim )
1533{
1534  pointSet *vs;
1535  onePoint vert;
1536  int j,k,l;
1537
1538  vert.point=(Coord_t*)Alloc( (pVariables+2) * sizeof(Coord_t) );
1539
1540  vs= new pointSet( dim );
1541
1542  for ( j= 1; j <= Q1->num; j++ )
1543  {
1544    for ( k= 1; k <= Q2->num; k++ )
1545    {
1546      for ( l= 1; l <= dim; l++ )
1547      {
1548        vert.point[l]= (*Q1)[j]->point[l] + (*Q2)[k]->point[l];
1549      }
1550      vs->mergeWithExp( &vert );
1551      //vs->addPoint( &vert );
1552    }
1553  }
1554
1555  Free( (ADDRESS) vert.point, (pVariables+2) * sizeof(Coord_t) );
1556
1557  return vs;
1558}
1559
1560pointSet * resMatrixSparse::minkSumAll( pointSet **pQ, int numq, int dim )
1561{
1562  pointSet *vs,*vs_old;
1563  int j;
1564
1565  vs= new pointSet( dim );
1566
1567  for ( j= 1; j <= pQ[0]->num; j++ ) vs->addPoint( (*pQ[0])[j] );
1568
1569  for ( j= 1; j < numq; j++ )
1570  {
1571    vs_old= vs;
1572    vs= minkSumTwo( vs_old, pQ[j], dim );
1573
1574    delete vs_old;
1575  }
1576
1577  return vs;
1578}
1579
1580//----------------------------------------------------------------------------------------
1581
1582resMatrixSparse::resMatrixSparse( const ideal _gls, const int special )
1583  : resMatrixBase(), gls( _gls )
1584{
1585  pointSet **Qi; // vertices sets of Conv(Supp(f_i)), i=0..idelem
1586  pointSet *E;   // all integer lattice points of the minkowski sum of Q0...Qn
1587  int i,j,k;
1588  int pnt;
1589  int totverts;                // total number of exponent vectors in ideal gls
1590  mprfloat shift[MAXVARS+2];   // shiftvector delta, index [1..dim]
1591
1592  if ( pVariables > MAXVARS )
1593  {
1594    WerrorS("resMatrixSparse::resMatrixSparse: To many variables!");
1595    return;
1596  }
1597
1598  rmat= NULL;
1599  numSet0= 0;
1600
1601  if ( special == SNONE ) linPolyS= 0;
1602  else linPolyS= special;
1603
1604  istate= resMatrixBase::ready;
1605
1606  n= pVariables;
1607  idelem= IDELEMS(gls);  // should be n+1
1608
1609  // prepare matrix LP->LiPM for Linear Programming
1610  totverts = 0;
1611  for( i=0; i < idelem; i++) totverts += pLength( (gls->m)[i] );
1612
1613  LP = new simplex( idelem+totverts*2+5, totverts+5 ); // rows, cols
1614
1615  // get shift vector
1616#ifdef mprTEST
1617  shift[0]=0.005; shift[1]=0.003; shift[2]=0.008; shift[3]=0.005; shift[4]=0.002;
1618  shift[5]=0.1; shift[6]=0.3; shift[7]=0.2; shift[8]=0.4; shift[9]=0.2;
1619#else
1620  randomVector( idelem, shift );
1621#endif
1622#ifdef mprDEBUG_PROT
1623  Print(" shift vector: ");
1624  for ( i= 1; i <= idelem; i++ ) Print(" %.12f ",(double)shift[i]);
1625  PrintLn();
1626#endif
1627
1628  // evaluate convex hull for supports of gls
1629  convexHull chnp( LP );
1630  Qi= chnp.newtonPolytopesP( gls );
1631
1632#ifdef mprMINKSUM
1633  E= minkSumAll( Qi, n+1, n);
1634#else
1635  // get inner points
1636  mayanPyramidAlg mpa( LP );
1637  E= mpa.getInnerPoints( Qi, shift );
1638#endif
1639
1640#ifdef mprDEBUG_PROT
1641#ifdef mprMINKSUM
1642  Print("(MinkSum)");
1643#endif
1644  PrintS("\n E = (Q_0 + ... + Q_n) \\cap \\N :\n");
1645  for ( pnt= 1; pnt <= E->num; pnt++ )
1646  {
1647    Print("%d: <",pnt);print_exp( (*E)[pnt], E->dim );Print(">\n");
1648  }
1649  PrintLn();
1650#endif
1651
1652#ifdef mprTEST
1653  int lift[5][5];
1654  lift[0][1]=3; lift[0][2]=4; lift[0][3]=8;  lift[0][4]=2;
1655  lift[1][1]=6; lift[1][2]=1; lift[1][3]=7;  lift[1][4]=4;
1656  lift[2][1]=2; lift[2][2]=5; lift[2][3]=9;  lift[2][4]=6;
1657  lift[3][1]=2; lift[3][2]=1; lift[3][3]=9;  lift[3][4]=5;
1658  lift[4][1]=3; lift[4][2]=7; lift[4][3]=1;  lift[4][4]=5;
1659  // now lift everything
1660  for ( i= 0; i <= n; i++ ) Qi[i]->lift( lift[i] );
1661  E->dim++;
1662#else
1663  // now lift everything
1664  for ( i= 0; i <= n; i++ ) Qi[i]->lift();
1665  E->dim++;
1666#endif
1667
1668  // run Row Content Function for every point in E
1669  for ( pnt= 1; pnt <= E->num; pnt++ )
1670  {
1671    RC( Qi, E, pnt, shift );
1672  }
1673
1674  // remove points not in cells
1675  k= E->num;
1676  for ( pnt= k; pnt > 0; pnt-- )
1677  {
1678    if ( (*E)[pnt]->rcPnt == NULL )
1679    {
1680      E->removePoint(pnt);
1681      mprSTICKYPROT(ST_SPARSE_RCRJ);
1682    }
1683  }
1684  mprSTICKYPROT("\n");
1685
1686#ifdef mprDEBUG_PROT
1687  PrintS(" points which lie in a cell:\n");
1688  for ( pnt= 1; pnt <= E->num; pnt++ )
1689  {
1690    Print("%d: <",pnt);print_exp( (*E)[pnt], E->dim );Print(">\n");
1691  }
1692  PrintLn();
1693#endif
1694
1695  // unlift to old dimension, sort
1696  for ( i= 0; i <= n; i++ ) Qi[i]->unlift();
1697  E->unlift();
1698  E->sort();
1699
1700#ifdef mprDEBUG_PROT
1701  Print(" points with a[ij] (%d):\n",E->num);
1702  for ( pnt= 1; pnt <= E->num; pnt++ )
1703  {
1704    Print("Punkt p \\in E[%d]: <",pnt);print_exp( (*E)[pnt], E->dim );
1705    Print(">, RC(p) = (i:%d, j:%d), a[i,j] = <",(*E)[pnt]->rc.set,(*E)[pnt]->rc.pnt);
1706    //print_exp( (Qi[(*E)[pnt]->rc.set])[(*E)[pnt]->rc.pnt], E->dim );Print("> = <");
1707    print_exp( (*E)[pnt]->rcPnt, E->dim );Print(">\n");
1708  }
1709#endif
1710
1711  // now create matrix
1712  if ( createMatrix( E ) != E->num )
1713  {
1714    // this can happen if the shiftvector shift is to large or not generic
1715    istate= resMatrixBase::fatalError;
1716    WerrorS("resMatrixSparse::resMatrixSparse: Error in resMatrixSparse::createMatrix!");
1717    goto theEnd;
1718    //return;
1719  }
1720
1721 theEnd:
1722  // clean up
1723  for ( i= 0; i < idelem; i++ )
1724  {
1725    delete Qi[i];
1726  }
1727  Free( (ADDRESS) Qi, idelem * sizeof(pointSet*) );
1728
1729  delete E;
1730
1731  delete LP;
1732}
1733
1734//----------------------------------------------------------------------------------------
1735
1736resMatrixSparse::~resMatrixSparse()
1737{
1738  delete uRPos;
1739  idDelete( &rmat );
1740}
1741
1742const ideal resMatrixSparse::getMatrix()
1743{
1744  int i,j,rp,cp;
1745  poly pp,phelp,piter,pgls;
1746
1747  // copy original sparse res matrix
1748  ideal rmat_out= idCopy(rmat);
1749
1750  // now fill in coeffs of f0
1751  for ( i= 1; i <= numSet0; i++ )
1752  {
1753
1754    pgls= (gls->m)[0]; // f0
1755
1756    // get matrix row and delete it
1757    pp= (rmat_out->m)[IMATELEM(*uRPos,i,1)];
1758    pDelete( &pp );
1759    pp= NULL;
1760    phelp= pp;
1761    piter= NULL;
1762
1763    // u_1,..,u_k
1764    cp=2;
1765    while ( pNext(pgls)!=NULL )
1766    {
1767      phelp= pOne();
1768      pSetCoeff( phelp, nCopy(pGetCoeff(pgls)) );
1769      pSetComp( phelp, IMATELEM(*uRPos,i,cp) );
1770      pSetmComp( phelp );
1771      if ( piter!=NULL )
1772      {
1773        pNext(piter)= phelp;
1774        piter= phelp;
1775      }
1776      else
1777      {
1778        pp= phelp;
1779        piter= phelp;
1780      }
1781      cp++;
1782      pIter( pgls );
1783    }
1784    // u0, now pgls points to last monom
1785    phelp= pOne();
1786    pSetCoeff( phelp, nCopy(pGetCoeff(pgls)) );
1787    //pSetComp( phelp, IMATELEM(*uRPos,i,idelem+1) );
1788    pSetComp( phelp, IMATELEM(*uRPos,i,pLength((gls->m)[0])+1) );
1789    pSetmComp( phelp );
1790    pNext(piter)= phelp;
1791    (rmat_out->m)[IMATELEM(*uRPos,i,1)]= pp;
1792  }
1793
1794  return rmat_out;
1795}
1796
1797// Fills in resMat[][] with evpoint[] and gets determinant
1798//    uRPos[i][1]: row of matrix
1799//    uRPos[i][idelem+1]: col of u(0)
1800//    uRPos[i][2..idelem]: col of u(1) .. u(n)
1801//    i= 1 .. numSet0
1802const number resMatrixSparse::getDetAt( const number* evpoint )
1803{
1804  int i,cp,rp;
1805  poly pp,phelp,piter;
1806
1807  mprPROTnl("smCallDet");
1808
1809  for ( i= 1; i <= numSet0; i++ )
1810  {
1811    pp= (rmat->m)[IMATELEM(*uRPos,i,1)];
1812    pDelete( &pp );
1813    pp= NULL;
1814    phelp= pp;
1815    piter= NULL;
1816    // u_1,..,u_n
1817    for ( cp= 2; cp <= idelem; cp++ )
1818    {
1819      if ( !nIsZero(evpoint[cp-1]) )
1820      {
1821        phelp= pOne();
1822        pSetCoeff( phelp, nCopy(evpoint[cp-1]) );
1823        pSetComp( phelp, IMATELEM(*uRPos,i,cp) );
1824        pSetmComp( phelp );
1825        if ( piter )
1826        {
1827          pNext(piter)= phelp;
1828          piter= phelp;
1829        }
1830        else
1831        {
1832          pp= phelp;
1833          piter= phelp;
1834        }
1835      }
1836    }
1837    // u0
1838    phelp= pOne();
1839    pSetCoeff( phelp, nCopy(evpoint[0]) );
1840    pSetComp( phelp, IMATELEM(*uRPos,i,idelem+1) );
1841    pSetmComp( phelp );
1842    pNext(piter)= phelp;
1843    (rmat->m)[IMATELEM(*uRPos,i,1)]= pp;
1844  }
1845
1846  mprSTICKYPROT(ST__DET); // 1
1847
1848  poly pres= smCallDet( rmat );
1849  number numres= nCopy( pGetCoeff( pres ) );
1850  pDelete( &pres );
1851
1852  mprSTICKYPROT(ST__DET); // 2
1853
1854  return ( numres );
1855}
1856
1857// Fills in resMat[][] with evpoint[] and gets determinant
1858//    uRPos[i][1]: row of matrix
1859//    uRPos[i][idelem+1]: col of u(0)
1860//    uRPos[i][2..idelem]: col of u(1) .. u(n)
1861//    i= 1 .. numSet0
1862const poly resMatrixSparse::getUDet( const number* evpoint )
1863{
1864  int i,cp,rp;
1865  poly pp,phelp,piter;
1866
1867  mprPROTnl("smCallDet");
1868
1869  for ( i= 1; i <= numSet0; i++ )
1870  {
1871    pp= (rmat->m)[IMATELEM(*uRPos,i,1)];
1872    pDelete( &pp );
1873    phelp= NULL;
1874    piter= NULL;
1875    for ( cp= 2; cp <= idelem; cp++ )
1876    { // u1 .. un
1877      if ( !nIsZero(evpoint[cp-1]) )
1878      {
1879        phelp= pOne();
1880        pSetCoeff( phelp, nCopy(evpoint[cp-1]) );
1881        pSetComp( phelp, IMATELEM(*uRPos,i,cp) );
1882        pSetmComp( phelp );
1883        if ( piter!=NULL )
1884        {
1885          pNext(piter)= phelp;
1886          piter= phelp;
1887        }
1888        else
1889        {
1890          pp= phelp;
1891          piter= phelp;
1892        }
1893      }
1894    }
1895    // u0
1896    phelp= pOne();
1897    pSetExp(phelp,1,1);
1898    pSetComp( phelp, IMATELEM(*uRPos,i,idelem+1) );
1899    pSetm( phelp );
1900    pNext(piter)= phelp;
1901    (rmat->m)[IMATELEM(*uRPos,i,1)]= pp;
1902  }
1903
1904  mprSTICKYPROT(ST__DET); // 1
1905
1906  poly pres= smCallDet( rmat );
1907
1908  mprSTICKYPROT(ST__DET); // 2
1909
1910  return ( pres );
1911}
1912//<-
1913
1914//-----------------------------------------------------------------------------
1915//-------------- dense resultant matrix ---------------------------------------
1916//-----------------------------------------------------------------------------
1917
1918//-> dense resultant matrix
1919//
1920class resVector;
1921
1922/* dense resultant matrix */
1923class resMatrixDense : virtual public resMatrixBase
1924{
1925public:
1926  /**
1927   * _gls: system of multivariate polynoms
1928   * special: -1 -> resMatrixDense is a symbolic matrix
1929   *    0,1, ... -> resMatrixDense ist eine u-Resultante, wobei special das
1930   *                        lineare u-Polynom angibt
1931   */
1932  resMatrixDense( const ideal _gls, const int special = SNONE );
1933  ~resMatrixDense();
1934
1935  /** column vector of matrix, index von 0 ... numVectors-1 */
1936  resVector *getMVector( const int i );
1937
1938  /** Returns the matrix M in an usable presentation */
1939  const ideal getMatrix();
1940
1941  /** Returns the submatrix M' of M in an usable presentation */
1942  const ideal getSubMatrix();
1943
1944  /** Evaluate the determinant of the matrix M at the point evpoint
1945   * where the ui's are replaced by the components of evpoint.
1946   * Uses singclap_det from fractory.
1947   */
1948  const number getDetAt( const number* evpoint );
1949
1950  /** Evaluates the determinant of the submatrix M'.
1951   * Since the matrix is numerically, no evaluation point is needed.
1952   * Uses singclap_det from fractory.
1953   */
1954  const number getSubDet();
1955
1956private:
1957  /** deactivated copy constructor */
1958  resMatrixDense( const resMatrixDense & );
1959
1960  /** Generate the "matrix" M. Each column is presented by a resVector
1961   * holding all entries for this column.
1962   */
1963  void generateBaseData();
1964
1965  /** Generates needed set of monoms, split them into sets S0, ... Sn and
1966   * check if reduced/nonreduced and calculate size of submatrix.
1967   */
1968  void generateMonomData( int deg, intvec* polyDegs , intvec* iVO );
1969
1970  /** Recursively generate all homogeneous monoms of
1971   * pVariables variables of degree deg.
1972   */
1973  void generateMonoms( poly m, int var, int deg );
1974
1975  /** Creates quadratic matrix M of size numVectors for later use.
1976   * u0, u1, ...,un are replaced by 0.
1977   * Entries equal to 0 are not initialized ( == NULL)
1978   */
1979  void createMatrix();
1980
1981private:
1982  resVector *resVectorList;
1983
1984  int veclistmax;
1985  int veclistblock;
1986  int numVectors;
1987  int subSize;
1988
1989  matrix m;
1990};
1991//<-
1992
1993//-> struct resVector
1994/* Holds a row vector of the dense resultant matrix */
1995struct resVector
1996{
1997public:
1998  void init()
1999  {
2000    isReduced = FALSE;
2001    elementOfS = SFREE;
2002    mon = NULL;
2003  }
2004  void init( const poly m )
2005  {
2006    isReduced = FALSE;
2007    elementOfS = SFREE;
2008    mon = m;
2009  }
2010
2011  /** index von 0 ... numVectors-1 */
2012  poly getElem( const int i );
2013
2014  /** index von 0 ... numVectors-1 */
2015  number getElemNum( const int i );
2016
2017  // variables
2018  poly mon;
2019  poly dividedBy;
2020  bool isReduced;
2021
2022  /** number of the set S mon is element of */
2023  int elementOfS;
2024
2025  /** holds the index of u0, u1, ..., un, if (elementOfS == linPolyS)
2026   *  the size is given by pVariables
2027   */
2028  int *numColParNr;
2029
2030  /** holds the column vector if (elementOfS != linPolyS) */
2031  number *numColVector;
2032
2033  /** size of numColVector */
2034  int numColVectorSize;
2035
2036  number *numColVecCopy;
2037};
2038//<-
2039
2040//-> resVector::*
2041poly resVector::getElem( const int i ) // inline ???
2042{
2043  assume( 0 < i || i > numColVectorSize );
2044  poly out= pOne();
2045  pSetCoeff( out, numColVector[i] );
2046  pTest( out );
2047  return( out );
2048}
2049
2050number resVector::getElemNum( const int i ) // inline ??
2051{
2052  assume( i >= 0 && i < numColVectorSize );
2053  return( numColVector[i] );
2054}
2055//<-
2056
2057//-> resMatrixDense::*
2058resMatrixDense::resMatrixDense( const ideal _gls, const int special )
2059  : resMatrixBase()
2060{
2061  int i;
2062
2063  sourceRing=currRing;
2064  gls= idCopy( _gls );
2065  linPolyS= special;
2066  m=NULL;
2067
2068  // init all
2069  generateBaseData();
2070
2071  totDeg= 1;
2072  for ( i= 0; i < IDELEMS(gls); i++ )
2073  {
2074    totDeg*=pTotaldegree( (gls->m)[i] );
2075  }
2076
2077  mprSTICKYPROT2("  resultant deg: %d\n",totDeg);
2078
2079  istate= resMatrixBase::ready;
2080}
2081
2082resMatrixDense::~resMatrixDense()
2083{
2084  int i,j;
2085  for (i=0; i < numVectors; i++)
2086  {
2087    pDelete( &resVectorList[i].mon );
2088    pDelete( &resVectorList[i].dividedBy );
2089    for ( j=0; j < resVectorList[i].numColVectorSize; j++ )
2090    {
2091        nDelete( resVectorList[i].numColVector+j );
2092    }
2093    Free( (ADDRESS)resVectorList[i].numColVector, numVectors*sizeof( number ) );
2094    Free( (ADDRESS)resVectorList[i].numColParNr, (pVariables+1) * sizeof(int) );
2095  }
2096
2097  Free( (ADDRESS)resVectorList, veclistmax*sizeof( resVector ) );
2098
2099  // free matrix m
2100  if ( m != NULL )
2101  {
2102    for ( i= 1; i <= numVectors; i++ )
2103      for ( j= 1; j <= numVectors; j++ )
2104        pDelete( &MATELEM(m , i, j) );
2105    Free( (ADDRESS)m->m, numVectors * numVectors * sizeof(poly) );
2106    FreeSizeOf( (ADDRESS)m, ip_smatrix );
2107  }
2108}
2109
2110// mprSTICKYPROT:
2111// ST_DENSE_FR: found row S0
2112// ST_DENSE_NR: normal row
2113void resMatrixDense::createMatrix()
2114{
2115  int k,i,j;
2116  resVector *vecp;
2117
2118  m= mpNew( numVectors, numVectors );
2119
2120  for ( i= 1; i <= MATROWS( m ); i++ )
2121    for ( j= 1; j <= MATCOLS( m ); j++ )
2122    {
2123      MATELEM(m,i,j)= pInit();
2124      pSetCoeff0( MATELEM(m,i,j), nInit(0) );
2125    }
2126
2127
2128  for ( k= 0; k <= numVectors - 1; k++ )
2129  {
2130    if ( linPolyS == getMVector(k)->elementOfS )
2131    {
2132      mprSTICKYPROT(ST_DENSE_FR);
2133      for ( i= 0; i < pVariables; i++ )
2134      {
2135        MATELEM(m,numVectors-k,numVectors-(getMVector(k)->numColParNr)[i])= pInit();
2136      }
2137    }
2138    else
2139    {
2140      mprSTICKYPROT(ST_DENSE_NR);
2141      vecp= getMVector(k);
2142      for ( i= 0; i < numVectors; i++)
2143      {
2144        if ( !nIsZero( vecp->getElemNum(i) ) )
2145        {
2146          MATELEM(m,numVectors - k,i + 1)= pInit();
2147          pSetCoeff0( MATELEM(m,numVectors - k,i + 1), nCopy(vecp->getElemNum(i)) );
2148        }
2149      }
2150    }
2151  } // for
2152  mprSTICKYPROT("\n");
2153
2154#ifdef mprDEBUG_ALL
2155  for ( k= numVectors - 1; k >= 0; k-- )
2156  {
2157    if ( linPolyS == getMVector(k)->elementOfS )
2158    {
2159      for ( i=0; i < pVariables; i++ )
2160      {
2161        Print(" %d ",(getMVector(k)->numColParNr)[i]);
2162      }
2163      PrintLn();
2164    }
2165  }
2166  for (i=1; i <= numVectors; i++)
2167  {
2168    for (j=1; j <= numVectors; j++ )
2169    {
2170      pWrite0(MATELEM(m,i,j));PrintS("  ");
2171    }
2172    PrintLn();
2173  }
2174#endif
2175}
2176
2177// mprSTICKYPROT:
2178// ST_DENSE_MEM: more mem allocated
2179// ST_DENSE_NMON: new monom added
2180void resMatrixDense::generateMonoms( poly m, int var, int deg )
2181{
2182  if ( deg == 0 )
2183  {
2184    poly mon = pCopy( m );
2185
2186    if ( numVectors == veclistmax )
2187    {
2188      resVectorList= (resVector * )ReAlloc( resVectorList,
2189                                            (veclistmax) * sizeof( resVector ),
2190                                            (veclistmax + veclistblock) * sizeof( resVector ) );
2191      int k;
2192      for ( k= veclistmax; k < (veclistmax + veclistblock); k++ )
2193        resVectorList[k].init();
2194      veclistmax+= veclistblock;
2195      mprSTICKYPROT(ST_DENSE_MEM);
2196
2197    }
2198    resVectorList[numVectors].init( mon );
2199    numVectors++;
2200    mprSTICKYPROT(ST_DENSE_NMON);
2201    return;
2202  }
2203  else
2204  {
2205    if ( var == pVariables+1 ) return;
2206    poly newm = pCopy( m );
2207    while ( deg >= 0 )
2208    {
2209      generateMonoms( newm, var+1, deg );
2210      pIncrExp( newm, var );
2211      pSetm( newm );
2212      deg--;
2213    }
2214    pDelete( & newm );
2215  }
2216
2217  return;
2218}
2219
2220void resMatrixDense::generateMonomData( int deg, intvec* polyDegs , intvec* iVO )
2221{
2222  int i,j,k;
2223
2224  // init monomData
2225  veclistblock= 512;
2226  veclistmax= veclistblock;
2227  resVectorList= (resVector *)Alloc( veclistmax*sizeof( resVector ) );
2228
2229  // Init resVector()s
2230  for ( j= veclistmax - 1; j >= 0; j-- ) resVectorList[j].init();
2231  numVectors= 0;
2232
2233  // Generate all monoms of degree deg
2234  poly start= pOne();
2235  generateMonoms( start, 1, deg );
2236  pDelete( & start );
2237
2238  mprSTICKYPROT("\n");
2239
2240  // Check for reduced monoms
2241  // First generate polyDegs.rows() monoms
2242  //  x(k)^(polyDegs[k]),  0 <= k < polyDegs.rows()
2243  ideal pDegDiv= idInit( polyDegs->rows(), 1 );
2244  for ( k= 0; k < polyDegs->rows(); k++ )
2245  {
2246    poly p= pOne();
2247    pSetExp( p, k + 1, (*polyDegs)[k] );
2248    pSetm( p );
2249    (pDegDiv->m)[k]= p;
2250  }
2251
2252  // Now check each monom if it is reduced.
2253  // A monom monom is called reduced if there exists
2254  // exactly one x(k)^(polyDegs[k]) that divides the monom.
2255  int divCount;
2256  for ( j= numVectors - 1; j >= 0; j-- )
2257  {
2258    divCount= 0;
2259    for ( k= 0; k < IDELEMS(pDegDiv); k++ )
2260      if ( pDivisibleBy2( (pDegDiv->m)[k], resVectorList[j].mon ) )
2261        divCount++;
2262    resVectorList[j].isReduced= (divCount == 1);
2263  }
2264
2265  // create the sets S(k)s
2266  // a monom x(i)^deg, deg given, is element of the set S(i)
2267  // if all x(0)^(polyDegs[0]) ... x(i-1)^(polyDegs[i-1]) DONT divide
2268  // x(i)^deg and only x(i)^(polyDegs[i]) divides x(i)^deg
2269  bool doInsert;
2270  for ( k= 0; k < iVO->rows(); k++)
2271  {
2272    //mprPROTInl(" ------------ var:",(*iVO)[k]);
2273    for ( j= numVectors - 1; j >= 0; j-- )
2274    {
2275      //mprPROTPnl("testing monom",resVectorList[j].mon);
2276      if ( resVectorList[j].elementOfS == SFREE )
2277      {
2278        //mprPROTnl("\tfree");
2279        if ( pDivisibleBy2( (pDegDiv->m)[ (*iVO)[k] ], resVectorList[j].mon ) )
2280        {
2281          //mprPROTPnl("\tdivisible by ",(pDegDiv->m)[ (*iVO)[k] ]);
2282          doInsert=TRUE;
2283          for ( i= 0; i < k; i++ )
2284          {
2285            //mprPROTPnl("\tchecking db ",(pDegDiv->m)[ (*iVO)[i] ]);
2286            if ( pDivisibleBy2( (pDegDiv->m)[ (*iVO)[i] ], resVectorList[j].mon ) )
2287            {
2288              //mprPROTPnl("\t and divisible by",(pDegDiv->m)[ (*iVO)[i] ]);
2289              doInsert=FALSE;
2290              break;
2291            }
2292          }
2293          if ( doInsert )
2294          {
2295            //mprPROTInl("\t------------------> S ",(*iVO)[k]);
2296            resVectorList[j].elementOfS= (*iVO)[k];
2297            resVectorList[j].dividedBy= pCopy( (pDegDiv->m)[ (*iVO)[i] ] );
2298          }
2299        }
2300      }
2301    }
2302  }
2303
2304  // size of submatrix M', equal to number of nonreduced monoms
2305  // (size of matrix M is equal to number of monoms=numVectors)
2306  subSize= 0;
2307  int sub;
2308  for ( i= 0; i < polyDegs->rows(); i++ )
2309  {
2310    sub= 1;
2311    for ( k= 0; k < polyDegs->rows(); k++ )
2312      if ( i != k ) sub*= (*polyDegs)[k];
2313    subSize+= sub;
2314  }
2315  subSize= numVectors - subSize;
2316
2317  // pDegDiv wieder freigeben!
2318  idDelete( &pDegDiv );
2319
2320#ifdef mprDEBUG_ALL
2321  // Print a list of monoms and their properties
2322  PrintS("// \n");
2323  for ( j= numVectors - 1; j >= 0; j-- )
2324  {
2325    Print("// %s, S(%d),  db ",
2326          resVectorList[j].isReduced?"reduced":"nonreduced",
2327          resVectorList[j].elementOfS);
2328    pWrite0(resVectorList[j].dividedBy);
2329    PrintS("  monom ");
2330    pWrite(resVectorList[j].mon);
2331  }
2332  Print("// size: %d, subSize: %d\n",numVectors,subSize);
2333#endif
2334}
2335
2336void resMatrixDense::generateBaseData()
2337{
2338  int k,j,i;
2339  number matEntry;
2340  poly pmatchPos;
2341  poly pi,factor,pmp;
2342
2343  // holds the degrees of F0, F1, ..., Fn
2344  intvec polyDegs( IDELEMS(gls) );
2345  for ( k= 0; k < IDELEMS(gls); k++ )
2346    polyDegs[k]= pTotaldegree( (gls->m)[k] );
2347
2348  // the internal Variable Ordering
2349  // make sure that the homogenization variable goes last!
2350  intvec iVO( pVariables );
2351  if ( linPolyS != SNONE )
2352  {
2353    iVO[pVariables - 1]= linPolyS;
2354    int p=0;
2355    for ( k= pVariables - 1; k >= 0; k-- )
2356    {
2357      if ( k != linPolyS )
2358      {
2359        iVO[p]= k;
2360        p++;
2361      }
2362    }
2363  }
2364  else
2365  {
2366    linPolyS= 0;
2367    for ( k= 0; k < pVariables; k++ )
2368      iVO[k]= pVariables - k - 1;
2369  }
2370
2371  // the critical degree d= sum( deg(Fi) ) - n
2372  int sumDeg= 0;
2373  for ( k= 0; k < polyDegs.rows(); k++ )
2374    sumDeg+= polyDegs[k];
2375  sumDeg-= polyDegs.rows() - 1;
2376
2377  // generate the base data
2378  generateMonomData( sumDeg, &polyDegs, &iVO );
2379
2380  // generate "matrix"
2381  for ( k= numVectors - 1; k >= 0; k-- )
2382  {
2383    if ( resVectorList[k].elementOfS != linPolyS )
2384    {
2385      // column k is a normal column with numerical or symbolic entries
2386      // init stuff
2387      resVectorList[k].numColParNr= NULL;
2388      resVectorList[k].numColVectorSize= numVectors;
2389      resVectorList[k].numColVector= (number *)Alloc( numVectors*sizeof( number ) );
2390      for ( i= 0; i < numVectors; i++ ) resVectorList[k].numColVector[i]= nInit(0);
2391
2392      // compute row poly
2393      poly pi= pMult( pCopy( (gls->m)[ resVectorList[k].elementOfS ] ) , pCopy( resVectorList[k].mon ) );
2394      pi= pDivideM( pCopy( pi ), pCopy( resVectorList[k].dividedBy ) );
2395
2396      // fill in "matrix"
2397      while ( pi != NULL )
2398      {
2399        matEntry= nCopy(pGetCoeff(pi));
2400        pmatchPos= pHead0( pi );
2401        pSetCoeff0( pmatchPos, nInit(1) );
2402
2403        for ( i= 0; i < numVectors; i++)
2404          if ( pEqual( pmatchPos, resVectorList[i].mon ) )
2405            break;
2406
2407        resVectorList[k].numColVector[numVectors - i - 1] = nCopy(matEntry);
2408
2409        pDelete( &pmatchPos );
2410        nDelete( &matEntry );
2411
2412        pIter( pi );
2413      }
2414      pDelete( &pi );
2415    }
2416    else
2417    {
2418      // column is a special column, i.e. is generated by S0 and F0
2419      // safe only the positions of the ui's in the column
2420      //mprPROTInl(" setup of numColParNr ",k);
2421      resVectorList[k].numColVectorSize= 0;
2422      resVectorList[k].numColVector= NULL;
2423      resVectorList[k].numColParNr= (int *)Alloc0( (pVariables+1) * sizeof(int) );
2424
2425      pi= (gls->m)[ resVectorList[k].elementOfS ];
2426      factor= pDivideM( pCopy( resVectorList[k].mon ), pCopy( resVectorList[k].dividedBy ) );
2427
2428      j=0;
2429      while ( pi  != NULL )
2430      { // fill in "matrix"
2431        pmp= pMult( pCopy( factor ), pHead( pi ) );
2432        pTest( pmp );
2433
2434        for ( i= 0; i < numVectors; i++)
2435          if ( pEqual( pmp, resVectorList[i].mon ) )
2436            break;
2437
2438        resVectorList[k].numColParNr[j]= i;
2439        pDelete( &pmp );
2440        pIter( pi );
2441        j++;
2442      }
2443      pDelete( &pi );
2444      pDelete( &factor );
2445    }
2446  } // for ( k= numVectors - 1; k >= 0; k-- )
2447
2448  mprSTICKYPROT2(" size of matrix:    %d\n",numVectors);
2449  mprSTICKYPROT2(" size of submatrix: %d\n",subSize);
2450
2451  // create the matrix M
2452  createMatrix();
2453}
2454
2455resVector *resMatrixDense::getMVector(int i)
2456{
2457  assume( i >= 0 && i < numVectors );
2458  return &resVectorList[i];
2459}
2460
2461const ideal resMatrixDense::getMatrix()
2462{
2463  int k,i,j;
2464
2465  // copy matrix
2466  matrix resmat= mpNew(numVectors,numVectors);
2467  for (i=1; i <= numVectors; i++)
2468  {
2469    for (j=1; j <= numVectors; j++ )
2470    {
2471      if ( (MATELEM(m,i,j)!=NULL)
2472      && (!nIsZero(pGetCoeff(MATELEM(m,i,j)))) )
2473      {
2474        MATELEM(resmat,i,j)= pCopy( MATELEM(m,i,j) );
2475      }
2476    }
2477  }
2478  for (i=0; i < numVectors; i++)
2479  {
2480    if ( resVectorList[i].elementOfS == linPolyS )
2481    {
2482      for (j=1; j <= pVariables; j++ )
2483      {
2484        if ( MATELEM(resmat,numVectors-i,
2485                     numVectors-resVectorList[i].numColParNr[j-1])!=NULL )
2486          pDelete( &MATELEM(resmat,numVectors-i,numVectors-resVectorList[i].numColParNr[j-1]) );
2487        MATELEM(resmat,numVectors-i,numVectors-resVectorList[i].numColParNr[j-1])= pOne();
2488        // FIX ME
2489        if ( true )
2490        {
2491          pSetCoeff( MATELEM(resmat,numVectors-i,numVectors-resVectorList[i].numColParNr[j-1]), nPar(j) );
2492        }
2493        else
2494        {
2495          pSetExp( MATELEM(resmat,numVectors-i,numVectors-resVectorList[i].numColParNr[j-1]), j, 1 );
2496          pSetm(MATELEM(resmat,numVectors-i,numVectors-resVectorList[i].numColParNr[j-1]));
2497        }
2498      }
2499    }
2500  }
2501
2502  // obachman: idMatrix2Module frees resmat !!
2503  ideal resmod= idMatrix2Module(resmat);
2504  return resmod;
2505}
2506
2507const ideal resMatrixDense::getSubMatrix()
2508{
2509  int k,i,j,l;
2510  resVector *vecp;
2511
2512  // generate quadratic matrix resmat of size subSize
2513  matrix resmat= mpNew( subSize, subSize );
2514
2515  j=1;
2516  for ( k= numVectors - 1; k >= 0; k-- )
2517  {
2518    vecp= getMVector(k);
2519    if ( vecp->isReduced ) continue;
2520    l=1;
2521    for ( i= numVectors - 1; i >= 0; i-- )
2522    {
2523      if ( getMVector(i)->isReduced ) continue;
2524      if ( !nIsZero(vecp->getElemNum(numVectors - i - 1)) )
2525      {
2526        MATELEM(resmat,j,l)= pCopy( vecp->getElem(numVectors-i-1) );
2527      }
2528      l++;
2529    }
2530    j++;
2531  }
2532
2533  // obachman: idMatrix2Module frees resmat !!
2534  ideal resmod= idMatrix2Module(resmat);
2535  return resmod;
2536}
2537
2538const number resMatrixDense::getDetAt( const number* evpoint )
2539{
2540  int k,i,j;
2541
2542  // copy evaluation point into matrix
2543  // p0, p1, ..., pn replace u0, u1, ..., un
2544  for ( k= numVectors - 1; k >= 0; k-- )
2545  {
2546    if ( linPolyS == getMVector(k)->elementOfS )
2547    {
2548      for ( i= 0; i < pVariables; i++ )
2549      {
2550        pSetCoeff( MATELEM(m,numVectors-k,numVectors-(getMVector(k)->numColParNr)[i]),
2551                   nCopy(evpoint[i]) );
2552      }
2553    }
2554  }
2555
2556  mprSTICKYPROT(ST__DET);
2557
2558  // evaluate determinant of matrix m using factory singclap_det
2559#ifdef HAVE_FACTORY
2560  poly res= singclap_det( m );
2561#else
2562  poly res= NULL;
2563#endif
2564
2565  // avoid errors for det==0
2566  number numres;
2567  if ( (res!=NULL)  && (!nIsZero(pGetCoeff( res ))) )
2568  {
2569    numres= nCopy( pGetCoeff( res ) );
2570  }
2571  else
2572  {
2573    numres= nInit(0);
2574    mprPROT("0");
2575  }
2576  pDelete( &res );
2577
2578  mprSTICKYPROT(ST__DET);
2579
2580  return( numres );
2581}
2582
2583const number resMatrixDense::getSubDet()
2584{
2585  int k,i,j,l;
2586  resVector *vecp;
2587
2588  // generate quadratic matrix mat of size subSize
2589  matrix mat= mpNew( subSize, subSize );
2590
2591  for ( i= 1; i <= MATROWS( mat ); i++ )
2592  {
2593    for ( j= 1; j <= MATCOLS( mat ); j++ )
2594    {
2595      MATELEM(mat,i,j)= pInit();
2596      pSetCoeff0( MATELEM(mat,i,j), nInit(0) );
2597    }
2598  }
2599  j=1;
2600  for ( k= numVectors - 1; k >= 0; k-- )
2601  {
2602    vecp= getMVector(k);
2603    if ( vecp->isReduced ) continue;
2604    l=1;
2605    for ( i= numVectors - 1; i >= 0; i-- )
2606    {
2607      if ( getMVector(i)->isReduced ) continue;
2608      if ( vecp->getElemNum(numVectors - i - 1) && !nIsZero(vecp->getElemNum(numVectors - i - 1)) )
2609      {
2610        pSetCoeff(MATELEM(mat, j , l ), nCopy(vecp->getElemNum(numVectors - i - 1)));
2611      }
2612      /* else
2613      {
2614           MATELEM(mat, j , l )= pOne();
2615           pSetCoeff(MATELEM(mat, j , l ), nInit(0) );
2616      }
2617      */
2618      l++;
2619    }
2620    j++;
2621  }
2622
2623#ifdef HAVE_FACTORY
2624  poly res= singclap_det( mat );
2625#else
2626  poly res= NULL;
2627#endif
2628
2629  number numres;
2630  if ((res != NULL) && (!nIsZero(pGetCoeff( res ))) )
2631  {
2632    numres= nCopy(pGetCoeff( res ));
2633  }
2634  else
2635  {
2636    numres= nInit(0);
2637  }
2638  pDelete( &res );
2639  return numres;
2640}
2641//<--
2642
2643//-----------------------------------------------------------------------------
2644//-------------- uResultant ---------------------------------------------------
2645//-----------------------------------------------------------------------------
2646
2647#define MAXEVPOINT 1000000 // 0x7fffffff
2648//#define MPR_MASI
2649
2650//-> unsigned long over(unsigned long n,unsigned long d)
2651// Calculates (n+d \over d) using gmp functionality
2652//
2653unsigned long over( const unsigned long n , const unsigned long d )
2654{ // (d+n)! / ( d! n! )
2655  mpz_t res;
2656  mpz_init(res);
2657  mpz_t m,md,mn;
2658  mpz_init(m);mpz_set_ui(m,1);
2659  mpz_init(md);mpz_set_ui(md,1);
2660  mpz_init(mn);mpz_set_ui(mn,1);
2661
2662  mpz_fac_ui(m,n+d);
2663  mpz_fac_ui(md,d);
2664  mpz_fac_ui(mn,n);
2665
2666  mpz_mul(res,md,mn);
2667  mpz_tdiv_q(res,m,res);
2668
2669  mpz_clear(m);mpz_clear(md);mpz_clear(mn);
2670
2671  unsigned long result = mpz_get_ui(res);
2672  mpz_clear(res);
2673
2674  return result;
2675}
2676//<-
2677
2678//-> uResultant::*
2679uResultant::uResultant( const ideal _gls, const resMatType _rmt, BOOLEAN extIdeal )
2680  : rmt( _rmt )
2681{
2682  if ( extIdeal )
2683  {
2684    // extend given ideal by linear poly F0=u0x0 + u1x1 +...+ unxn
2685    gls= extendIdeal( _gls, linearPoly( rmt ), rmt );
2686    n= IDELEMS( gls );
2687  }
2688  else
2689    gls= idCopy( _gls );
2690
2691  switch ( rmt )
2692  {
2693  case sparseResMat:
2694    resMat= new resMatrixSparse( gls );
2695    break;
2696  case denseResMat:
2697#ifdef HAVE_FACTORY
2698    resMat= new resMatrixDense( gls );
2699    break;
2700#endif
2701  default:
2702    WerrorS("uResultant::uResultant: Unknown resultant matrix type choosen!");
2703  }
2704}
2705
2706uResultant::~uResultant( )
2707{
2708  delete resMat;
2709}
2710
2711ideal uResultant::extendIdeal( const ideal gls, poly linPoly, const resMatType rmt )
2712{
2713  ideal newGls= idCopy( gls );
2714  newGls->m= (poly *)ReAlloc( newGls->m,
2715                              IDELEMS(gls) * sizeof(poly),
2716                              (IDELEMS(gls) + 1) * sizeof(poly) );
2717  IDELEMS(newGls)++;
2718
2719  switch ( rmt )
2720  {
2721  case sparseResMat:
2722  case denseResMat:
2723    {
2724      int i;
2725      for ( i= IDELEMS(newGls)-1; i > 0; i-- )
2726      {
2727        newGls->m[i]= newGls->m[i-1];
2728      }
2729      newGls->m[0]= linPoly;
2730    }
2731    break;
2732  default:
2733    WerrorS("uResultant::extendIdeal: Unknown resultant matrix type choosen!");
2734  }
2735
2736  return( newGls );
2737}
2738
2739poly uResultant::linearPoly( const resMatType rmt )
2740{
2741  int i,j;
2742
2743  poly newlp= pOne();
2744  poly actlp, rootlp= newlp;
2745
2746  for ( i= 1; i <= pVariables; i++ )
2747  {
2748    actlp= newlp;
2749    pSetExp( actlp, i, 1 );
2750    pSetm( actlp );
2751    newlp= pOne();
2752    actlp->next= newlp;
2753  }
2754  actlp->next= NULL;
2755  pDelete( &newlp );
2756
2757  if ( rmt == sparseResMat )
2758  {
2759    newlp= pOne();
2760    actlp->next= newlp;
2761    newlp->next= NULL;
2762  }
2763  return ( rootlp );
2764}
2765
2766poly uResultant::interpolateDense( const number subDetVal )
2767{
2768  int i,j,p;
2769  long tdg;
2770
2771  // D is a Polynom homogeneous in the coeffs of F0 and of degree tdg = d0*d1*...*dn
2772  tdg= resMat->getDetDeg();
2773
2774  // maximum number of terms in polynom D (homogeneous, of degree tdg)
2775  // long mdg= (facul(tdg+n-1) / facul( tdg )) / facul( n - 1 );
2776  long mdg= over( n-1, tdg );
2777
2778  // maxiaml number of terms in a polynom of degree tdg
2779  long l=(long)pow( tdg+1, n );
2780
2781#ifdef mprDEBUG_PROT
2782  Print("// total deg of D: tdg %d\n",tdg);
2783  Print("// maximum number of terms in D: mdg: %d\n",mdg);
2784  Print("// maximum number of terms in polynom of deg tdg: l %d\n",l);
2785#endif
2786
2787  // we need mdg results of D(p0,p1,...,pn)
2788  number *presults;
2789  presults= (number *)Alloc( mdg * sizeof( number ) );
2790  for (i=0; i < mdg; i++) presults[i]= nInit(0);
2791
2792  number *pevpoint= (number *)Alloc( n * sizeof( number ) );
2793  number *pev= (number *)Alloc( n * sizeof( number ) );
2794  for (i=0; i < n; i++) pev[i]= nInit(0);
2795
2796  mprPROTnl("// initial evaluation point: ");
2797  // initial evaluatoin point
2798  p=1;
2799  for (i=0; i < n; i++)
2800  {
2801    // init pevpoint with primes 3,5,7,11, ...
2802    p= nextPrime( p );
2803    pevpoint[i]=nInit( p );
2804    nTest(pevpoint[i]);
2805    mprPROTNnl(" ",pevpoint[i]);
2806  }
2807
2808  // evaluate the determinant in the points pev^0, pev^1, ..., pev^mdg
2809  mprPROTnl("// evaluating:");
2810  for ( i=0; i < mdg; i++ )
2811  {
2812    for (j=0; j < n; j++)
2813    {
2814      nDelete( &pev[j] );
2815      nPower(pevpoint[j],i,&pev[j]);
2816      mprPROTN(" ",pev[j]);
2817    }
2818    mprPROTnl("");
2819
2820    nDelete( &presults[i] );
2821    presults[i]=resMat->getDetAt( pev );
2822
2823    mprSTICKYPROT(ST_BASE_EV);
2824  }
2825  mprSTICKYPROT("\n");
2826
2827  // now interpolate using vandermode interpolation
2828  mprPROTnl("// interpolating:");
2829  number *ncpoly;
2830  {
2831    vandermonde vm( mdg, n, tdg, pevpoint );
2832    ncpoly= vm.interpolateDense( presults );
2833  }
2834
2835  if ( subDetVal != NULL )
2836  {   // divide by common factor
2837    number detdiv;
2838    for ( i= 0; i <= mdg; i++ )
2839    {
2840      detdiv= nDiv( ncpoly[i], subDetVal );
2841      nNormalize( detdiv );
2842      nDelete( &ncpoly[i] );
2843      ncpoly[i]= detdiv;
2844    }
2845  }
2846
2847#ifdef mprDEBUG_ALL
2848  PrintLn();
2849  for ( i=0; i < mdg; i++ )
2850  {
2851    nPrint(ncpoly[i]); PrintS(" --- ");
2852  }
2853  PrintLn();
2854#endif
2855
2856  // prepare ncpoly for later use
2857  number nn=nInit(0);
2858  for ( i=0; i < mdg; i++ )
2859  {
2860    if ( nEqual(ncpoly[i],nn) )
2861    {
2862      nDelete( &ncpoly[i] );
2863      ncpoly[i]=NULL;
2864    }
2865  }
2866  nDelete( &nn );
2867
2868  // create poly presenting the determinat of the uResultant
2869  intvec exp( n );
2870  for ( i= 0; i < n; i++ ) exp[i]=0;
2871
2872  poly result= NULL;
2873
2874  long sum=0;
2875  long c=0;
2876
2877  for ( i=0; i < l; i++ )
2878  {
2879    if ( sum == tdg )
2880    {
2881      if ( !nIsZero(ncpoly[c]) )
2882      {
2883        poly p= pOne();
2884        if ( rmt == denseResMat )
2885        {
2886          for ( j= 0; j < n; j++ ) pSetExp( p, j+1, exp[j] );
2887        }
2888        else if ( rmt == sparseResMat )
2889        {
2890          for ( j= 1; j < n; j++ ) pSetExp( p, j, exp[j] );
2891        }
2892        pSetCoeff( p, ncpoly[c] );
2893        pSetm( p );
2894        if (result!=NULL) result= pAdd( result, p );
2895        else result=  p;
2896      }
2897      c++;
2898    }
2899    sum=0;
2900    exp[0]++;
2901    for ( j= 0; j < n - 1; j++ )
2902    {
2903      if ( exp[j] > tdg )
2904      {
2905        exp[j]= 0;
2906        exp[j + 1]++;
2907      }
2908      sum+=exp[j];
2909    }
2910    sum+=exp[n-1];
2911  }
2912
2913  pTest( result );
2914
2915  return result;
2916}
2917
2918rootContainer ** uResultant::interpolateDenseSP( BOOLEAN matchUp, const number subDetVal )
2919{
2920  int i,j,p,uvar;
2921  long tdg;
2922  int loops= (matchUp?n-2:n-1);
2923
2924  mprPROTnl("uResultant::interpolateDenseSP");
2925
2926  tdg= resMat->getDetDeg();
2927
2928  // evaluate D in tdg+1 distinct points, so
2929  // we need tdg+1 results of D(p0,1,0,...,0) =
2930  //              c(0)*u0^tdg + c(1)*u0^tdg-1 + ... + c(tdg-1)*u0 + c(tdg)
2931  number *presults;
2932  presults= (number *)Alloc( (tdg + 1) * sizeof( number ) );
2933  for ( i=0; i <= tdg; i++ ) presults[i]= nInit(0);
2934
2935  rootContainer ** roots;
2936  roots= (rootContainer **) Alloc( loops * sizeof(rootContainer*) );
2937  for ( i=0; i < loops; i++ ) roots[i]= new rootContainer(); // 0..n-2
2938
2939  number *pevpoint= (number *)Alloc( n * sizeof( number ) );
2940  for (i=0; i < n; i++) pevpoint[i]= nInit(0);
2941
2942  number *pev= (number *)Alloc( n * sizeof( number ) );
2943  for (i=0; i < n; i++) pev[i]= nInit(0);
2944
2945  // now we evaluate D(u0,-1,0,...0), D(u0,0,-1,0,...,0), ..., D(u0,0,..,0,-1)
2946  // or D(u0,k1,k2,0,...,0), D(u0,k1,k2,k3,0,...,0), ..., D(u0,k1,k2,k3,...,kn)
2947  // this gives us n-1 evaluations
2948  p=3;
2949  for ( uvar= 0; uvar < loops; uvar++ )
2950  {
2951    // generate initial evaluation point
2952    if ( matchUp )
2953    {
2954      for (i=0; i < n; i++)
2955      {
2956        // prime(random number) between 1 and MAXEVPOINT
2957        nDelete( &pevpoint[i] );
2958        if ( i == 0 )
2959        {
2960          //p= nextPrime( p );
2961          pevpoint[i]= nInit( p );
2962        }
2963        else if ( i <= uvar + 2 )
2964        {
2965          pevpoint[i]=nInit(1+siRand()%MAXEVPOINT);
2966          //pevpoint[i]=nInit(383);
2967        }
2968        else
2969          pevpoint[i]=nInit(0);
2970        mprPROTNnl(" ",pevpoint[i]);
2971      }
2972    }
2973    else
2974    {
2975      for (i=0; i < n; i++)
2976      {
2977        // init pevpoint with  prime,0,...0,1,0,...,0
2978        nDelete( &pevpoint[i] );
2979        if ( i == 0 )
2980        {
2981          //p=nextPrime( p );
2982          pevpoint[i]=nInit( p );
2983        }
2984        else
2985        {
2986          if ( i == (uvar + 1) ) pevpoint[i]= nInit(-1);
2987          else pevpoint[i]= nInit(0);
2988        }
2989        mprPROTNnl(" ",pevpoint[i]);
2990      }
2991    }
2992
2993    // prepare aktual evaluation point
2994    for (i=0; i < n; i++)
2995    {
2996      nDelete( &pev[i] );
2997      pev[i]= nCopy( pevpoint[i] );
2998    }
2999    // evaluate the determinant in the points pev^0, pev^1, ..., pev^tdg
3000    for ( i=0; i <= tdg; i++ )
3001    {
3002      nDelete( &pev[0] );
3003      nPower(pevpoint[0],i,&pev[0]);          // new evpoint
3004
3005      nDelete( &presults[i] );
3006      presults[i]=resMat->getDetAt( pev );   // evaluate det at point evpoint
3007
3008      mprPROTNnl("",presults[i]);
3009
3010      mprSTICKYPROT(ST_BASE_EV);
3011      mprPROTI("",tdg-i);
3012    }
3013    mprSTICKYPROT("\n");
3014
3015    // now interpolate
3016    vandermonde vm( tdg + 1, 1, tdg, pevpoint, FALSE );
3017    number *ncpoly= vm.interpolateDense( presults );
3018
3019    if ( subDetVal != NULL )
3020    {  // divide by common factor
3021      number detdiv;
3022      for ( i= 0; i <= tdg; i++ )
3023      {
3024        detdiv= nDiv( ncpoly[i], subDetVal );
3025        nNormalize( detdiv );
3026        nDelete( &ncpoly[i] );
3027        ncpoly[i]= detdiv;
3028      }
3029    }
3030
3031#ifdef mprDEBUG_ALL
3032    PrintLn();
3033    for ( i=0; i <= tdg; i++ )
3034    {
3035      nPrint(ncpoly[i]); PrintS(" --- ");
3036    }
3037    PrintLn();
3038#endif
3039
3040    // save results
3041    roots[uvar]->fillContainer( ncpoly, pevpoint, uvar+1, tdg,
3042                                (matchUp?rootContainer::cspecialmu:rootContainer::cspecial),
3043                                loops );
3044  }
3045
3046  // free some stuff: pev, presult
3047  for ( i=0; i < n; i++ ) nDelete( pev + i );
3048  Free( (ADDRESS)pev, n * sizeof( number ) );
3049
3050  for ( i=0; i <= tdg; i++ ) nDelete( presults+i );
3051  Free( (ADDRESS)presults, (tdg + 1) * sizeof( number ) );
3052
3053  return roots;
3054}
3055
3056rootContainer ** uResultant::specializeInU( BOOLEAN matchUp, const number subDetVal )
3057{
3058  int i,j,p,uvar;
3059  long tdg;
3060  poly pures,piter;
3061  int loops=(matchUp?n-2:n-1);
3062
3063  mprPROTnl("uResultant::specializeInU");
3064
3065  tdg= resMat->getDetDeg();
3066
3067  rootContainer ** roots;
3068  roots= (rootContainer **) Alloc( loops * sizeof(rootContainer*) );
3069  for ( i=0; i < loops; i++ ) roots[i]= new rootContainer(); // 0..n-2
3070
3071  number *pevpoint= (number *)Alloc( n * sizeof( number ) );
3072  for (i=0; i < n; i++) pevpoint[i]= nInit(0);
3073
3074  // now we evaluate D(u0,-1,0,...0), D(u0,0,-1,0,...,0), ..., D(u0,0,..,0,-1)
3075  // or D(u0,k1,k2,0,...,0), D(u0,k1,k2,k3,0,...,0), ..., D(u0,k1,k2,k3,...,kn)
3076  p=3;
3077  for ( uvar= 0; uvar < loops; uvar++ )
3078  {
3079    // generate initial evaluation point
3080    if ( matchUp )
3081    {
3082      for (i=0; i < n; i++)
3083      {
3084        // prime(random number) between 1 and MAXEVPOINT
3085        nDelete( &pevpoint[i] );
3086        if ( i <= uvar + 2 )
3087        {
3088          pevpoint[i]=nInit(1+siRand()%MAXEVPOINT);
3089          //pevpoint[i]=nInit(383);
3090        } else pevpoint[i]=nInit(0);
3091        mprPROTNnl(" ",pevpoint[i]);
3092      }
3093    }
3094    else
3095    {
3096      for (i=0; i < n; i++)
3097      {
3098        // init pevpoint with  prime,0,...0,-1,0,...,0
3099        nDelete( &(pevpoint[i]) );
3100        if ( i == (uvar + 1) ) pevpoint[i]= nInit(-1);
3101        else pevpoint[i]= nInit(0);
3102        mprPROTNnl(" ",pevpoint[i]);
3103      }
3104    }
3105
3106    pures= resMat->getUDet( pevpoint );
3107
3108    number *ncpoly= (number *)Alloc( (tdg+1) * sizeof(number) );
3109
3110#ifdef MPR_MASI
3111    BOOLEAN masi=true;
3112#endif
3113
3114    piter= pures;
3115    for ( i= tdg; i >= 0; i-- )
3116    {
3117      //if ( piter ) Print("deg %d, pDeg(piter) %d\n",i,pTotaldegree(piter));
3118      if ( piter && pTotaldegree(piter) == i )
3119      {
3120        ncpoly[i]= nCopy( pGetCoeff( piter ) );
3121        pIter( piter );
3122#ifdef MPR_MASI
3123        masi=false;
3124#endif
3125      }
3126      else
3127      {
3128        ncpoly[i]= nInit(0);
3129      }
3130      mprPROTNnl("", ncpoly[i] );
3131    }
3132#ifdef MPR_MASI
3133    if ( masi ) mprSTICKYPROT("MASI");
3134#endif
3135
3136    mprSTICKYPROT(ST_BASE_EV); // .
3137
3138    if ( subDetVal != NULL ) {  // divide by common factor
3139      number detdiv;
3140      for ( i= 0; i <= tdg; i++ )
3141      {
3142        detdiv= nDiv( ncpoly[i], subDetVal );
3143        nNormalize( detdiv );
3144        nDelete( &ncpoly[i] );
3145        ncpoly[i]= detdiv;
3146      }
3147    }
3148
3149    pDelete( &pures );
3150
3151    // save results
3152    roots[uvar]->fillContainer( ncpoly, pevpoint, uvar+1, tdg,
3153                                (matchUp?rootContainer::cspecialmu:rootContainer::cspecial),
3154                                loops );
3155  }
3156
3157  mprSTICKYPROT("\n");
3158
3159  // free some stuff: pev, presult
3160  for ( i=0; i < n; i++ ) nDelete( pevpoint + i );
3161  Free( (ADDRESS)pevpoint, n * sizeof( number ) );
3162
3163  return roots;
3164}
3165
3166int uResultant::nextPrime( int i )
3167{
3168  int init=i;
3169  i+=2;
3170  int j= IsPrime( i );
3171  while ( j <= init )
3172  {
3173    i+=2;
3174    j= IsPrime( i );
3175  }
3176  return j;
3177}
3178//<-
3179
3180//-----------------------------------------------------------------------------
3181
3182//-> loNewtonPolytope(...)
3183ideal loNewtonPolytope( const ideal id )
3184{
3185  simplex * LP;
3186  int i;
3187  int n,totverts,idelem;
3188  ideal idr;
3189
3190  n= pVariables;
3191  idelem= IDELEMS(id);  // should be n+1
3192
3193  totverts = 0;
3194  for( i=0; i < idelem; i++) totverts += pLength( (id->m)[i] );
3195
3196  LP = new simplex( idelem+totverts*2+5, totverts+5 ); // rows, cols
3197
3198  // evaluate convex hull for supports of id
3199  convexHull chnp( LP );
3200  idr = chnp.newtonPolytopesI( id );
3201
3202  delete LP;
3203
3204  return idr;
3205}
3206//<-
3207
3208//%e
3209
3210//-----------------------------------------------------------------------------
3211
3212#endif // ! ASO_GENERATE
3213
3214// local Variables: ***
3215// folded-file: t ***
3216// compile-command-1: "make installg" ***
3217// compile-command-2: "make install" ***
3218// End: ***
3219
3220// in folding: C-c x
3221// leave fold: C-c y
3222//   foldmode: F10
Note: See TracBrowser for help on using the repository browser.