source: git/coeffs/mpr_base.cc @ 655a29

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