source: git/Singular/mpr_base.cc @ 35aab3

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