source: git/Singular/mpr_base.cc @ 2a7c6d

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