source: git/numeric/mpr_base.cc @ a9c298

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