Gröbner basis project
Codebase for research into Gröbner basis computation
algorithm_buchberger_explorer.cpp
1 #ifndef __ALGORITHM_BUCHBERGER_EXPLORER_CPP_
2 #define __ALGORITHM_BUCHBERGER_EXPLORER_CPP_
3 
4 /*****************************************************************************\
5 * This file is part of DynGB. *
6 * *
7 * DynGB is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 2 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * Foobar is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have received a copy of the GNU General Public License *
18 * along with DynGB. If not, see <http://www.gnu.org/licenses/>. *
19 \*****************************************************************************/
20 
21 #include <vector>
22 using std::vector;
23 
24 #include "hilbert_functions.hpp"
25 
26 #include "algorithm_buchberger_basic.hpp"
27 
28 #include "reduction_support.hpp"
29 
30 #include "algorithm_buchberger_explorer.hpp"
31 
32 #include "mpi.h"
33 
34 #define XPLOR_GENERATOR -1
35 #define PROC_UNASSIGNED -1
36 
37 #define HILBERT_COMPARISON 0
38 #define SUGAR_COMPARISON 1
39 #define WHICH_COMPARISON SUGAR_COMPARISON
40 
47 public:
49 
52  : Critical_Pair_Basic(f, strategy)
53  {
54  pi = i; qi = XPLOR_GENERATOR; proc = PROC_UNASSIGNED;
55  }
58  int i, int j, StrategyFlags strategy, vector<Abstract_Polynomial *>G
59  )
60  : Critical_Pair_Basic(G[i], G[j], strategy)
61  {
62  pi = i; qi = j; proc = PROC_UNASSIGNED;
63  }
66  int i, Abstract_Polynomial *g, StrategyFlags strategy,
67  vector<Abstract_Polynomial *>G
68  )
69  : Critical_Pair_Basic(G[i], g, strategy)
70  {
71  pi = i; qi = G.size(); proc = PROC_UNASSIGNED;
72  }
74 
75 
77  int first_index() { return pi; }
79  int second_index() { return qi; }
81  DEG_TYPE sugar() {
82  return (static_cast<Pair_Sugar_Data *>(key))->pair_sugar();
83  }
85 
86 
88  void set_processor(int i) { proc = i; }
93  int get_processor() { return proc; }
95 protected:
97  int pi;
99  int qi;
101  int proc;
102 };
103 
110 bool second_HF_smaller(
115 ) {
116  // this adapts LessByHilbert in dynamic_engine.cpp
117  // in this case, we want the opposite result,
118  // so the return statement at the end is negated
119  bool result;
120  // first check the coefficients of the Hilbert polynomial
122  HPdiff -= *hp2;
123  if (not HPdiff.is_zero())
124  result = (HPdiff.numerator(HPdiff.degree()) < 0);
125  else // use Hilbert series
126  {
129  DEG_TYPE i = 0;
130  for ( /* already initialized */ ;
131  i <= h1->degree() and i <= h2->degree() and (*h1)[i] == (*h2)[i];
132  i++)
133  { /* taken care of in loop */ }
134  if (i > h1->degree())
135  {
136  if (i > h2->degree())
137  // the numerators are equal; second is not measurably better
138  result = false;
139  else
140  result = true;
141  }
142  else
143  {
144  if (i > h2->degree()) result = false;
145  else result = (*h1)[i] < (*h2)[i];
146  }
147  }
148  //cout << "\tfirst less than second? " << result << endl;
149  return not result;
150 }
151 
158 bool newer_HF_smaller(Monomial & t, unsigned i, Monomial & u, unsigned j,
161 ) {
162  // this adapts LessByHilbert in dynamic_engine.cpp
163  // in this case, we want the opposite result,
164  // so the return statement at the end is negated
165  bool result;
166  NVAR_TYPE n = t.num_vars();
167  // first check the coefficients of the Hilbert polynomial
168  Dense_Univariate_Rational_Polynomial HPdiff(*(HP[i]));
169  HPdiff -= *(HP[j]);
170  if (not HPdiff.is_zero())
171  result = (HPdiff.numerator(HPdiff.degree()) < 0);
172  else // use Hilbert series
173  {
176  DEG_TYPE i = 0;
177  for ( /* already initialized */ ;
178  i <= h1->degree() and i <= h2->degree() and (*h1)[i] == (*h2)[i];
179  i++)
180  { /* taken care of in loop */ }
181  if (i > h1->degree())
182  {
183  if (i > h2->degree())
184  { // the numerators are equal; break tie via lex
185  int i = 0;
186  while (i < n and t[i] == u[i]) ++i;
187  if (i == n) result = false;
188  else result = (t.degree(i) > u.degree(i));
189  }
190  else
191  result = true;
192  }
193  else
194  {
195  if (i > h2->degree()) result = false;
196  else result = (*h1)[i] < (*h2)[i];
197  }
198  }
199  //cout << "\tfirst less than second? " << result << endl;
200  return not result;
201 }
202 
217 void hilbert_chooser(
218  void * firstin, void * secondin, int * len, MPI_Datatype * dptr
219 ) {
220  int64_t * first = (int64_t *)firstin;
221  int64_t * second = (int64_t *)secondin;
222  bool undecided = true;
223  bool keep_second = true;
224  if (first[0] == -1 or second[0] == -1) {
225  undecided = false;
226  if (first[0] != -1)
227  keep_second = false;
228  }
229  if (undecided) {
230  int64_t i;
231  int64_t n = first[1];
232  for (i = 0; undecided and i <= n; ++i)
233  undecided = (first[2*i + 2]*second[2*i+3] == second[2*i + 2]*first[2*i+3]);
234  if (not undecided) {
235  --i;
236  keep_second = (first[2*i + 2]*second[2*i+3] >= second[2*i + 2]*first[2*i+3]);
237  }
238  else {
239  for (i = 2*(n + 1) + 2; undecided and i < *len; ++i)
240  undecided = (first[i] == second[i]);
241  --i;
242  keep_second = undecided or first[i] < second[i];
243  }
244  }
245  // at this point, if we're still undecided then the functions are equal
246  if (not keep_second)
247  for (unsigned i = 0; i < *len; ++i)
248  second[i] = first[i];
249 }
250 
255 void sugar_chooser(
256  void * firstin, void * secondin, int * len, MPI_Datatype * dptr
257 ) {
258  int64_t * first = (int64_t *)firstin;
259  int64_t * second = (int64_t *)secondin;
260  int64_t n = first[2];
261  bool keep_second = true;
262  if (first[1] == -1) {
263  /* do nothing */
264  } else if (second[1] == -1 or first[2] > second[2]) {
265  keep_second = false;
266  } else if (first[1] > second[1]) {
267  /* do nothing */
268  // this comparison MUST be separated from the above
269  // otherwise (nonnegative, -1) makes us keep -1
270  } else {
271  Monomial t(n);
272  Monomial u(n);
273  for (unsigned i = 0; i < n; ++i) {
274  t.set_exponent(i, first[i+3]);
275  u.set_exponent(i, second[i+3]);
276  }
277  keep_second = t > u;
278  }
279  if (not keep_second) {
280  for (unsigned i = 0; i < n + 3; ++i)
281  second[i] = first[i];
282  }
283 }
284 
297  list<Critical_Pair_XPlor *> & P,
298  list<Critical_Pair_XPlor *> & Pass,
299  list<Critical_Pair_XPlor *> * Pcancel,
300  vector<Abstract_Polynomial *> & G,
302  StrategyFlags strategy
303 ) {
304  //cout << "----------------------\n";
305  list<Critical_Pair_XPlor *> C;
306  //cout << "creating and pruning pairs, starting with:\n";
307  //cout << '\t' << P.size() << " pairs\n";
308  //cout << '\t' << G.size() << " polynomials\n";
309  // critical pairs with new polynomial
310  unsigned m = G.size();
311  for (unsigned i = 0; i < G.size(); ++i)
312  C.push_back(new Critical_Pair_XPlor(i, r, strategy, G));
313  //cout << "Created " << C.size() << " pairs with new polynomial\n";
314  // apply Buchberger's lcm criterion to new pairs
315  list<Critical_Pair_XPlor *> D;
316  while (C.size() != 0) {
317  Critical_Pair_XPlor * p = C.front();
318  C.pop_front();
319  if ((p->first()->leading_monomial().is_coprime(
320  p->second()->leading_monomial()))
321  or (no_triplet(p, C) and no_triplet(p, D))
322  )
323  D.push_back(p);
324  else {
325  //cout << "triplet prunes " << *p << endl;
326  delete p;
327  }
328  }
329  //cout << "After applying Buchberger's lcm criterion to new pairs, now have "
330  // << D.size() << " new pairs\n";
331  // apply Buchberger's gcd criterion
332  list<Critical_Pair_XPlor *> E;
333  while (D.size() != 0) {
334  Critical_Pair_XPlor * p = D.front();
335  D.pop_front();
336  if (!(p->first()->leading_monomial().is_coprime(
337  p->second()->leading_monomial())))
338  E.push_back(p);
339  else {
340  //cout << "gcd prunes " << *p << endl;
341  delete p;
342  }
343  }
344  //cout << "After applying Buchberger's gcd criterion to new pairs, now have "
345  // << E.size() << " new pairs\n";
346  // apply Buchberger's lcm criterion to old pairs
347  list<Critical_Pair_XPlor *> Q;
348  while (P.size() != 0) {
349  Critical_Pair_XPlor * p = P.front();
350  P.pop_front();
351  if (!(r->leading_monomial() | p->lcm())
354  )
355  Q.push_back(p);
356  else {
357  //cout << "triplet prunes " << *p << endl;
358  delete p;
359  }
360  }
361  //cout << "After applying Buchberger's lcm criterion to new pairs, now have "
362  // << Q.size() << " old pairs\n";
363  P = Q;
364  list <Critical_Pair_XPlor *> R;
365  while (Pass.size() != 0) {
366  Critical_Pair_XPlor * p = Pass.front();
367  Pass.pop_front();
368  if (!(r->leading_monomial() | p->lcm())
371  )
372  R.push_back(p);
373  else {
374  //cout << "\ttriplet prunes " << *p << endl;
375  Pcancel[p->get_processor()].push_back(p);
376  }
377  }
378  Pass = R;
379  // add new pairs to old pairs
380  for (Critical_Pair_XPlor * e : E)
381  P.push_back(e);
382  /*for (unsigned i = 0; i < 2; ++i) {
383  if (Pcancel[i].size() > 0) {
384  cout << "Should delete pairs ";
385  for (Critical_Pair_XPlor * p : Pcancel[i])
386  cout << '(' << p->first_index() << ',' << p->second_index() << " (" << i << ") ) ";
387  cout << endl;
388  }
389  }*/
390  /*cout << "All pairs:\n";
391  for (list<Critical_Pair_XPlor *>::iterator pi = P.begin(); pi != P.end(); ++pi)
392  cout << '\t' << **pi << endl;
393  cout << "----------------------\n";*/
394 }
395 
397 
401 typedef struct {
403  int first;
405  int second;
407  DEG_TYPE sugar;
409 
410 list<Constant_Polynomial *> buchberger_explorer(
411  const vector<Abstract_Polynomial *> &F,
412  SPolyCreationFlags method,
413  StrategyFlags strategy,
414  WT_TYPE * strategy_weights,
415  const int comm_id,
416  const int comm_size
417 ) {
418  double r_bcast_time = 0;
419  unsigned number_of_spolys = 0;
420  vector<Abstract_Polynomial *> G; // basis
421  list<Critical_Pair_XPlor *> P; // critical pairs
422  list<Critical_Pair_XPlor *> Pass, * Pdel; // critical pair assignments
423  list<Critical_Pair_XPlor *> R; // reduced polynomials
424  list<Constant_Polynomial *> B; // end result
425  Polynomial_Ring & Rx = (*(F.begin()))->base_ring();
426  const Monomial_Ordering * mord = (*(F.begin()))->monomial_ordering();
427  // set up MPI_Datatype
428  MPI_Datatype pair_type, pair_basetypes[2];
429  MPI_Aint pair_offsets[2], extent, lb;
430  int pair_block_counts[2];
431  pair_offsets[0] = 0; pair_basetypes[0] = MPI_INT; pair_block_counts[0] = 2;
432  MPI_Type_get_extent(MPI_INT, &lb, &extent);
433  pair_offsets[1] = 2 * extent; pair_basetypes[1] = MPI_UNSIGNED_LONG_LONG;
434  pair_block_counts[1] = 1;
435  MPI_Type_create_struct(2, pair_block_counts, pair_offsets, pair_basetypes, &pair_type);
436  MPI_Type_commit(&pair_type);
437  // set up MPI_Op
438  MPI_Op comparison_op;
439  #if (WHICH_COMPARISON == HILBERT_COMPARISON)
440  MPI_Op_create(hilbert_chooser, false, &comparison_op);
441  #else
442  MPI_Op_create(sugar_chooser, false, &comparison_op);
443  #endif
444  // set up basis with generators
445  vector<Abstract_Polynomial *>::const_iterator Fi = F.begin();
446  NVAR_TYPE n = (*Fi)->number_of_variables();
447  for (unsigned i = 0; i < F.size(); ++i, ++Fi)
448  {
449  Abstract_Polynomial * fo = *Fi;
451  f->set_strategy(new Poly_Sugar_Data(f));
452  if (f->strategy() != nullptr) { f->strategy()->at_generation_tasks(); }
453  if (comm_id == 0) // only control needs to take care of critical pairs
454  P.push_back(new Critical_Pair_XPlor(i, strategy, f));
455  }
456  // main loop
457  bool verbose = false;
458  bool very_verbose = false;
459  unsigned min_todo = 0;
460  if (comm_id == 0) {
461  Pdel = new list<Critical_Pair_XPlor *>[comm_size];
463  = new Dense_Univariate_Rational_Polynomial *[comm_size];
465  = new Dense_Univariate_Integer_Polynomial *[comm_size];
467  = new Dense_Univariate_Integer_Polynomial *[comm_size];
468  min_todo = P.size();
469  }
470  MPI_Bcast(&min_todo, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
471  list<Monomial> T;
473  // create, reduce s-polynomials
476  int number_in_package;
477  while (min_todo != 0) {
478  /*for (int i = 0; i < comm_size; ++i) {
479  if (comm_id == i) {
480  cout << comm_id << "'s pairs\n";
481  for (Critical_Pair_XPlor * p : R)
482  cout << '\t' << comm_id << ' ' << *p << endl;
483  }
484  MPI_Barrier(MPI_COMM_WORLD);
485  }*/
486  //double start_time = MPI_Wtime();
487  uint64_t my_sugar = (R.size() == 0) ? 0 : R.front()->sugar();
488  uint64_t curr_sugar;
489  MPI_Reduce(&my_sugar, &curr_sugar, 1, MPI_UINT64_T, MPI_MAX, 0, MPI_COMM_WORLD);
490  if (comm_id == 0) {
491  if (P.size() > 1)
493  cout << "estimate " << min_todo << " steps left\n";
494  // first select comm_id pairs for reduction, send to coprocessors, reduce
495  int i = 1;
497  if (curr_sugar == 0 and !P.empty()) {
498  curr_sugar = P.front()->sugar();
499  for (list<Critical_Pair_XPlor *>::iterator pi = P.begin(); pi != P.end(); ++pi)
500  if ((*pi)->sugar() < curr_sugar)
501  curr_sugar = (*pi)->sugar();
502  }
503  cout << "current sugar: " << curr_sugar << endl;
504  unsigned num_of_curr_sugar = 0;
505  for (list<Critical_Pair_XPlor *>::iterator pi = P.begin();
506  pi != P.end();
507  ++pi)
508  {
509  if ((*pi)->sugar() == curr_sugar) ++num_of_curr_sugar;
510  }
511  number_in_package = 0;
512  if (num_of_curr_sugar == 0)
513  for (unsigned i = 1; i < comm_size; ++i)
514  MPI_Send(&number_in_package, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
515  else {
516  unsigned distrib_size = num_of_curr_sugar / comm_size + 1;
517  Pnew = new Critical_Pair_Communication[distrib_size];
518  list<Critical_Pair_XPlor *> Q;
519  for (unsigned i = 1; i < comm_size; ++i) {
520  number_in_package = 0;
521  while (!P.empty() and number_in_package < distrib_size)
522  {
523  Critical_Pair_XPlor * p = P.front();
524  P.pop_front();
525  if (p->sugar() != curr_sugar)
526  Q.push_back(p);
527  else {
528  Pnew[number_in_package].first = p->first_index();
529  Pnew[number_in_package].second = p->second_index();
530  Pnew[number_in_package].sugar = p->sugar();
531  ++number_in_package;
532  p->set_processor(i);
533  Pass.push_back(p);
534  }
535  }
536  MPI_Send(&number_in_package, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
537  if (number_in_package > 0)
538  MPI_Send(Pnew, number_in_package, pair_type, i, 0, MPI_COMM_WORLD);
539  }
540  number_in_package = 0;
541  while (!P.empty() and number_in_package < distrib_size)
542  {
543  Critical_Pair_XPlor * p = P.front();
544  P.pop_front();
545  if (p->sugar() != curr_sugar)
546  Q.push_back(p);
547  else {
548  Pnew[number_in_package].first = p->first_index();
549  Pnew[number_in_package].second = p->second_index();
550  Pnew[number_in_package].sugar = p->sugar();
551  ++number_in_package;
552  p->set_processor(0);
553  Pass.push_back(p);
554  }
555  }
556  P = Q;
557  }
558  /*for ( ; !P.empty() and i < comm_size; ++i) {
559  Critical_Pair_XPlor * p = P.front();
560  p_new.first = p->first_index(); p_new.second = p->second_index();
561  p_new.sugar = p->sugar();
562  //report_front_pair(p, strategy);
563  P.pop_front();
564  p->set_processor(i);
565  Pass.push_back(p);
566  //cout << comm_id << " sending "
567  // << p_new.first << ',' << p_new.second << ':' << p_new.sugar
568  // << " to " << i << endl;
569  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
570  }*/
571  /*for ( ; P.empty() and i < comm_size; ++i) {
572  p_new.first = p_new.second = XPLOR_GENERATOR; p_new.sugar = 0;
573  cout << comm_id << " sending "
574  // << p_new.first << ',' << p_new.second << ':' << p_new.sugar
575  // << " to " << i << endl;
576  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
577  }*/
578  /*if (!P.empty()) {
579  Critical_Pair_XPlor * p = P.front();
580  P.pop_front();
581  //report_front_pair(p, strategy);
582  p_in.first = p->first_index(); p_in.second = p->second_index();
583  p_in.sugar = p->sugar();
584  p->set_processor(comm_id);
585  Pass.push_back(p);
586  //cout << comm_id << " sending "
587  // << p_in.first << ',' << p_in.second << ':' << p_in.sugar
588  // << " to " << comm_id << endl;
589  } else {
590  p_in.first = p_in.second = XPLOR_GENERATOR;
591  } */
592  }
593  MPI_Barrier(MPI_COMM_WORLD);
594  if (comm_id != 0) {
595  MPI_Recv(&number_in_package, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
596  if (number_in_package > 0) {
597  Pnew = new Critical_Pair_Communication[number_in_package];
598  MPI_Recv(Pnew, number_in_package, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
599  }
600  }
601  for (unsigned i = 0; i < number_in_package; ++i) {
602  Critical_Pair_Communication & Pin = Pnew[i];
603  if (Pin.second == XPLOR_GENERATOR)
604  R.push_back(new Critical_Pair_XPlor(
605  Pin.first, StrategyFlags::SUGAR_STRATEGY, F[Pin.first]
606  ));
607  else
608  R.push_back(new Critical_Pair_XPlor(
609  Pin.first, Pin.second, StrategyFlags::SUGAR_STRATEGY, G
610  ));
611  }
612  if (number_in_package > 0)
613  delete [] Pnew;
614  /*if (comm_id != 0)
615  MPI_Recv(&p_in, 1, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
616  if (p_in.first != XPLOR_GENERATOR) {
617  if (p_in.second == XPLOR_GENERATOR)
618  //R.push_back(new Critical_Pair_XPlor(p_in.first, p_in.second, p_in.sugar, F));
619  R.push_back(new Critical_Pair_XPlor(p_in.first, StrategyFlags::SUGAR_STRATEGY, F[p_in.first]));
620  else
621  R.push_back(new Critical_Pair_XPlor(p_in.first, p_in.second, StrategyFlags::SUGAR_STRATEGY, G));
622  } */
623  /*for (unsigned i = 0; i < comm_size; ++i) {
624  if (comm_id == i) {
625  for (list<Critical_Pair_XPlor *>::iterator Ri = R.begin(); Ri != R.end(); ++Ri) {
626  cout << comm_id << " has " << (*Ri)->first_index() << ',' << (*Ri)->second_index() << " : " << (*Ri)->sugar() << endl;
627  }
628  }
629  MPI_Barrier(MPI_COMM_WORLD);
630  }
631  if (comm_id == 0) cout << "----\n";*/
632  //r_bcast_time += MPI_Wtime() - start_time;
633  //MPI_Barrier(MPI_COMM_WORLD);
634  double start_time = MPI_Wtime();
635  for (Critical_Pair_XPlor * p : R) {
636  // make s-poly
637  if ((s = p->s_polynomial()) == nullptr) {
638  s = p->s_polynomial(method, strategy);
639  ++number_of_spolys;
640  }
641  //cout << comm_id << " reducing s-poly "
642  // << p->first_index() << ',' << p->second_index() << endl;
643  //double start_time = MPI_Wtime();
644  //if (!s->is_zero() and find_reducer(s, G))
645  // reduce_over_basis<vector<Abstract_Polynomial *>>(&s, G, comm_id);
647  while (!s->is_zero() and ((g = find_reducer(s, G)) != nullptr)) {
648  //cout << comm_id << " reducing s-poly " << *s << endl;
649  top_reduce(s, g);
650  //cout << comm_id << " reduced to " << *s << endl;
651  }
652  //r_bcast_time += MPI_Wtime() - start_time;
653  p->set_spoly(s);
654  }
655  r_bcast_time += MPI_Wtime() - start_time;
656  /*for (int i = 0; i < comm_size; ++i) {
657  if (comm_id == i)
658  cout << comm_id << " finished reduction with " << R.size() << " polynomials\n";
659  MPI_Barrier(MPI_COMM_WORLD);
660  }*/
661  // create, compare Hilbert functions
662  //MPI_Barrier(MPI_COMM_WORLD);
663  //double start_time = MPI_Wtime();
664  unsigned winning_index = R.size() + 1;
665  Monomial * wt = nullptr;
666  #if (WHICH_COMPARISON == HILBERT_COMPARISON)
667  Dense_Univariate_Rational_Polynomial * wHP = nullptr;
668  Dense_Univariate_Integer_Polynomial * wHS = nullptr;
670  Dense_Univariate_Integer_Polynomial * hsn = nullptr;
672  #else
673  int64_t wd = -1;
674  #endif
675  // loop through critical pairs to find optimal pair
676  list<Critical_Pair_XPlor *>::iterator Rbest = R.end();
677  //double start_time = MPI_Wtime();
678  for (
679  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
680  Ri != R.end();
681  /* advance manually */
682  ) {
683  s = (*Ri)->s_polynomial();
684  if (s->is_zero()) {
685  //cout << '\t' << comm_id << ' ' << (*Ri)->first_index() << ','
686  // << (*Ri)->second_index() << ": reduced to zero\n";
687  Critical_Pair_XPlor * p = *Ri;
688  if (Ri == R.begin()) {
689  R.pop_front();
690  Ri = R.begin();
691  } else {
692  list<Critical_Pair_XPlor *>::iterator Rdel = Ri;
693  ++Ri;
694  R.erase(Rdel);
695  }
696  delete p;
697  } else {
698  #if (WHICH_COMPARISON == HILBERT_COMPARISON)
699  T.push_back(s->leading_monomial());
700  unsigned n = T.front().num_vars();
702  hsn = hilbert_second_numerator(n, hn);
703  unsigned d = ideal_dimension(n, hn, hsn);
704  hp = hilbert_polynomial(n, d, T, hn, hsn);
705  T.pop_back();
706  if (wt == nullptr) {
707  wt = &(s->leading_monomial()); wHP = hp; wHS = hn; Rbest = Ri;
708  } else {
709  if (second_HF_smaller(wHP, wHS, hp, hn)) {
710  wt = &(s->leading_monomial()); wHP = hp; wHS = hn; Rbest = Ri;
711  } else {
712  //delete hp;
713  //delete hn;
714  }
715  }
716  #else
717  if (wt == nullptr) {
718  wt = &(s->leading_monomial()); Rbest = Ri;
719  wd = (static_cast<Poly_Sugar_Data *>(s->strategy()))->poly_sugar();
720  } else {
721  if (*wt > s->leading_monomial()) {
722  wt = &(s->leading_monomial()); Rbest = Ri;
723  wd = (static_cast<Poly_Sugar_Data *>(s->strategy()))->poly_sugar();
724  }
725  }
726  #endif
727  ++Ri;
728  }
729  }
730  //r_bcast_time = MPI_Wtime() - start_time;
731  //MPI_Barrier(MPI_COMM_WORLD);
732  // move result of winning reduction to basis; return others to P
733  //cout << comm_id << " about to send Hilbert data\n";
734  /*for (int i = 0; i < comm_size; ++i) {
735  if (comm_id == i) {
736  if (wt == nullptr) {
737  cout << comm_id << " has no Hilbert data\n";
738  } else {
739  cout << comm_id << "'s Hilbert data for " << *wt << ":\n";
740  cout << '\t' << *wHP << endl;
741  cout << '\t' << *wHS << endl;
742  }
743  }
744  MPI_Barrier(MPI_COMM_WORLD);
745  }*/
746  /*for (int i = 0; i < comm_size; ++i) {
747  if (comm_id == i)
748  cout << comm_id << " sending " << wd << endl;
749  MPI_Barrier(MPI_COMM_WORLD);
750  }*/
751  //MPI_Barrier(MPI_COMM_WORLD);
752  //double start_time = MPI_Wtime();
753  int winners_id = -1;
754  // we use -1 if there is no useful data; i.e., all R's polys -> 0
755  #if (WHICH_COMPARISON == HILBERT_COMPARISON)
756  int hp_deg = (wt == nullptr) ? -1 : hp->degree();
757  int recv_hp_deg;
758  MPI_Allreduce(&hp_deg, &recv_hp_deg, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
759  // it's possible all the polynomials reduced to zero, so...
760  if (recv_hp_deg != -1) {
761  int hn_deg = (wHP == nullptr) ? -1 : wHS->degree();
762  int recv_hn_deg;
763  MPI_Allreduce(&hn_deg, &recv_hn_deg, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
764  int xfer_size = 2*(recv_hp_deg +1) + (recv_hn_deg + 1) + 2;
765  int64_t * send_data = new int64_t[xfer_size];
766  if (wt == nullptr)
767  send_data[0] = -1;
768  else {
769  send_data[0] = comm_id; send_data[1] = recv_hp_deg;
770  DEG_TYPE hi = 0;
771  for ( /* already initialized */ ; hi < recv_hp_deg - wHP->degree(); ++hi)
772  send_data[2*hi + 2] = send_data[2*hi + 3] = 0;
773  for ( /* already initialized */ ; hi <= recv_hp_deg; ++hi) {
774  send_data[2*hi + 2] = wHP->numerator(recv_hp_deg - hi);
775  send_data[2*hi + 3] = wHP->denominator(recv_hp_deg - hi);
776  }
777  for (hi = 0; hi <= wHS->degree(); ++hi)
778  send_data[2 + 2*(recv_hp_deg + 1) + hi] = wHS->coeff(hi);
779  for (
780  hi = 2 + 2*(recv_hp_deg + 1) + wHS->degree() + 1;
781  hi < xfer_size;
782  ++hi
783  ) {
784  send_data[hi] = 0;
785  }
786  }
787  int64_t * recv_data = new int64_t[xfer_size];
788  MPI_Allreduce(
789  send_data, recv_data, xfer_size, MPI_INT64_T, comparison_op,
790  MPI_COMM_WORLD
791  );
792  winners_id = recv_data[0];
793  delete [] send_data; delete [] recv_data;
794  }
795  #else
796  int64_t * send_data = new int64_t[n + 3];
797  int64_t * recv_data = new int64_t[n + 3];
798  send_data[0] = comm_id; send_data[1] = wd; send_data[2] = n;
799  if (wd >= 0) {
800  for (unsigned i = 0; i < n; ++i)
801  send_data[i + 3] = wt->degree(i);
802  }
803  MPI_Allreduce(
804  send_data, recv_data, n + 3, MPI_INT64_T, comparison_op,
805  MPI_COMM_WORLD
806  );
807  winners_id = (recv_data[1] == -1) ? -1 : recv_data[0];
808  delete [] send_data; delete [] recv_data;
809  //cout << "winner of comparison: " << winners_id << endl;
810  #endif
811  //delete wHP; delete wHS;
812  //r_bcast_time += MPI_Wtime() - start_time;
813  //if (comm_id == 0) cout << "chose winner " << winners_id << endl;
814  //double start_time = MPI_Wtime();
815  if (winners_id >= 0) {
816  uint64_t * r_bcast;
817  uint64_t r_bcast_size;
818  uint64_t r_bcast_sugar;
820  if (comm_id == winners_id) {
821  // notify proc 0
822  int pair_data[2];
823  pair_data[0] = (*Rbest)->first_index();
824  pair_data[1] = (*Rbest)->second_index();
825  if (comm_id != 0)
826  MPI_Send(pair_data, 2, MPI_INT, 0, 0, MPI_COMM_WORLD);
827  // broadcast winning polynomial
828  s = (*Rbest)->s_polynomial();
829  reduce_over_basis(&s, G, comm_id);
830  r = new Constant_Polynomial(*s);
831  Poly_Sugar_Data * sd = static_cast<Poly_Sugar_Data *>(s->strategy());
832  r->set_strategy(sd);
833  cout << "selected " << (*Rbest)->first_index()
834  << ',' << (*Rbest)->second_index() << ": " << r->leading_monomial()
835  << ',' << sd->poly_sugar() << " from " << winners_id << endl;
836  r_bcast_sugar = sd->poly_sugar();
837  s->set_strategy(nullptr);
838  delete s;
839  Critical_Pair_XPlor * p = *Rbest;
840  R.erase(Rbest);
841  delete p;
842  r_bcast = r->serialized(r_bcast_size);
843  }
844  if (comm_id == 0 and winners_id != 0) {
845  int pair_data[2];
846  MPI_Recv(pair_data, 2, MPI_INT, winners_id, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
847  bool found = false;
848  list<Critical_Pair_XPlor *>::iterator Ri;
849  for (Ri = Pass.begin(); not found and Ri != Pass.end(); /* manual */) {
850  if ((*Ri)->first_index() == pair_data[0] and (*Ri)->second_index() == pair_data[1])
851  found = true;
852  else
853  ++Ri;
854  }
855  if (found) {
856  //cout << "0 deleting " << pair_data[0] << ',' << pair_data[1] << endl;
857  delete *Ri;
858  Pass.erase(Ri);
859  }
860  }
861  //MPI_Barrier(MPI_COMM_WORLD);
862  MPI_Bcast(&r_bcast_size, 1, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
863  if (comm_id != winners_id)
864  r_bcast = new uint64_t [r_bcast_size];
865  //double start_time = MPI_Wtime();
866  MPI_Bcast(r_bcast, r_bcast_size, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
867  //r_bcast_time += MPI_Wtime() - start_time;
868  MPI_Bcast(&r_bcast_sugar, 1, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
869  // other processors need to create a copy of the new polynomial
870  if (comm_id != winners_id) {
871  r_bcast_size /= n + 1; // adjust size from # of words to # of terms
872  r = new Constant_Polynomial(Rx, mord, r_bcast_size, r_bcast);
873  Poly_Sugar_Data * rd = new Poly_Sugar_Data(r);
874  rd->force_sugar(r_bcast_sugar);
875  r->set_strategy(rd);
876  }
877  //if (comm_id == 0) cout << "selected " << r->leading_monomial() << endl;
878  //MPI_Barrier(MPI_COMM_WORLD);
879  //r_bcast_time += MPI_Wtime() - start_time;
880  delete [] r_bcast;
881  //double start_time = MPI_Wtime();
882  if (comm_id == 0) {
883  //cout << "added " << G.size() << ": " << r->leading_monomial() << endl;
884  very_verbose = false;
885  if (very_verbose) { cout << "\tin full "; r->println(); }
886  very_verbose = false;
887  gm_update_explorer(P, Pass, Pdel, G, r, strategy);
888  // remove useless pairs that were sent out to processors
889  for (int i = 1; i < comm_size; ++i) {
890  int outgoing = Pdel[i].size();
891  MPI_Send(&outgoing, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
892  while (Pdel[i].size() > 0) {
893  Critical_Pair_XPlor * p = Pdel[i].front();
895  Pdel[i].pop_front();
896  p_new.first = p->first_index();
897  p_new.second = p->second_index();
898  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
899  delete p;
900  }
901  }
902  // now delete my own redundant pairs
903  while (Pdel[0].size() > 0) {
904  Critical_Pair_XPlor * p = Pdel[0].front();
905  Pdel[0].pop_front();
906  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
907  while (
908  Ri != R.end() and
909  ((*Ri)->first_index() != p->first_index() or
910  (*Ri)->second_index() != p->second_index())
911  )
912  ++Ri;
913  // if the polynomial reduced to zero, it will not be in R
914  if (Ri != R.end()) {
915  //cout << comm_id << " deleting " << '(' << p->first_index() << ',' << p->second_index() << ')' << endl;
916  R.erase(Ri);
917  }
918  delete p;
919  }
920  } else {
921  int incoming;
922  MPI_Recv(&incoming, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
923  for (int j = 0; j < incoming; ++j) {
924  MPI_Recv(&p_in, 1, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
925  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
926  while (
927  Ri != R.end() and
928  ((*Ri)->first_index() != p_in.first or
929  (*Ri)->second_index() != p_in.second)
930  )
931  ++Ri;
932  // if Ri has reduced to zero, it will not be in R
933  if (Ri != R.end()) {
934  Critical_Pair_XPlor * p = *Ri;
935  //cout << comm_id << " deleting " << '(' << p->first_index() << ',' << p->second_index() << ')' << endl;
936  R.erase(Ri);
937  delete p;
938  }
939  }
940  }
941  //r_bcast_time += MPI_Wtime() - start_time;
942  //MPI_Barrier(MPI_COMM_WORLD);
943  G.push_back(r);
944  T.push_back(r->leading_monomial());
945  }
946  //r_bcast_time += MPI_Wtime() - start_time;
947  /*for (int i = 0; i < comm_size; ++i) {
948  if (comm_id == i) {
949  cout << comm_id << "'s basis\n";
950  for (Abstract_Polynomial * g : G)
951  cout << '\t' << comm_id << ' ' << *g << endl;
952  }
953  MPI_Barrier(MPI_COMM_WORLD);
954  }*/
955  //double start_time = MPI_Wtime();
956  if (comm_id == 0) min_todo = P.size();
957  unsigned my_todo = R.size();
958  unsigned all_todo;
959  MPI_Reduce(&my_todo, &all_todo, 1, MPI_UNSIGNED, MPI_MAX, 0, MPI_COMM_WORLD);
960  if (comm_id == 0)
961  min_todo = min_todo < all_todo ? all_todo : min_todo;
962  MPI_Bcast(&min_todo, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
963  //r_bcast_time += MPI_Wtime() - start_time;
964  //cout << comm_id << " understands there to be " << min_todo << " pairs remaining\n";
965  /*for (int i = 0; i < comm_size; ++i) {
966  if (comm_id == i)
967  for (Critical_Pair_XPlor * p : R)
968  cout << comm_id << "has " << p->first_index() << ',' << p->second_index() << endl;
969  MPI_Barrier(MPI_COMM_WORLD);
970  }*/
971  //MPI_Barrier(MPI_COMM_WORLD);
972  }
973  if (comm_id != 0) {
974  for (Abstract_Polynomial * g : G)
975  delete g;
976  }
977  MPI_Type_free(&pair_type);
978  cout << comm_id << " reduced " << number_of_spolys << endl;
979  //MPI_Barrier(MPI_COMM_WORLD);
980  int total_spolys;
981  MPI_Reduce(&number_of_spolys, &total_spolys, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
982  double max_bcast_time;
983  MPI_Reduce(&r_bcast_time, &max_bcast_time, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
984  for (unsigned i = 0; i < comm_size; ++i) {
985  if (comm_id == i)
986  cout << comm_id << " took " << r_bcast_time << " milliseconds on timed segment(s).\n";
987  MPI_Barrier(MPI_COMM_WORLD);
988  }
989  if (comm_id == 0) {
990  cout << total_spolys << " s-polynomials computed and reduced\n";
991  cout << max_bcast_time << " seconds spent on timed segment(s)\n";
992  // cleanup
993  cout << G.size() << " polynomials before interreduction\n";
994  //check_correctness(G, strategy);
995  list<Abstract_Polynomial *> G_final;
996  for (Abstract_Polynomial * g : G)
997  G_final.push_back(g);
998  G_final = reduce_basis(G_final);
999  cout << G_final.size() << " polynomials after interreduction\n";
1000  //set<Constant_Polynomial *, smaller_lm> B;
1001  for (Abstract_Polynomial * g : G_final) {
1002  B.push_back(new Constant_Polynomial(*g));
1003  //if (F.find(g) == F.end()) delete g;
1004  }
1005  delete [] Pdel;
1006  }
1007  return B;
1008 }
1009 
1010 #endif
The general class of a polynomial.
Definition: polynomial.hpp:101
list< Constant_Polynomial * > buchberger_explorer(const vector< Abstract_Polynomial *> &F, SPolyCreationFlags method, StrategyFlags strategy, WT_TYPE *strategy_weights, const int comm_id, const int comm_size)
Alternate implementation of Buchberger’s algorithm, for parallelization.
virtual bool is_zero() const =0
is this polynomial zero?
MPZCOEF_TYPE numerator(DEG_TYPE k) const
returns the th numerator
int second
index of second polynomial in the pair
A Constant_Polynomial is a polynomial that should not change.
void sort_pairs_by_strategy(list< T *> &P)
Applies the strategy to find the “smallest” critical pair.
int first
index of first polynomial in the pair
Dense_Univariate_Rational_Polynomial * hilbert_polynomial(NVAR_TYPE n, unsigned int pole_order, const list< Monomial > T, Dense_Univariate_Integer_Polynomial *hn, Dense_Univariate_Integer_Polynomial *hn2)
computes the Hilbert polynomial for an ideal
bool is_zero() const
indicates whether the polynomial is zero
Critical_Pair_XPlor(int i, int j, StrategyFlags strategy, vector< Abstract_Polynomial *>G)
create critical pair (f,g) where f, g are at indices i, j
const Abstract_Polynomial * second() const
second polynomial in pair
void gm_update_explorer(list< Critical_Pair_XPlor *> &P, list< Critical_Pair_XPlor *> &Pass, list< Critical_Pair_XPlor *> *Pcancel, vector< Abstract_Polynomial *> &G, Abstract_Polynomial *r, StrategyFlags strategy)
Implementation of Gebauer-Moeller algorithm, with XPLOR critical pairs. Based on description in Becke...
Mutable_Polynomial * s
S-polynomial.
int first_index()
returns index of first polynomial in pair
used to pass inforation on a critical pair from one polynomial to another
virtual Monomial & leading_monomial() const =0
leading monomial – call after sort_by_order()!
const Monomial & lcm() const
lcm of leading monomials of polynomials
int get_processor()
query whether this pair is assigned to a processor, and which (nonnegative value indicates assignment...
Dense_Univariate_Integer_Polynomial * hilbert_second_numerator(NVAR_TYPE n, Dense_Univariate_Integer_Polynomial *first, const WT_TYPE *grading)
computes the second Hilbert numerator (after reduction by )
void set_exponent(NVAR_TYPE i, DEG_TYPE e)
change th exponent to
Definition: monomial.cpp:78
void set_processor(int i)
record that this pair is assigned to processor i
int second_index()
returns index of second polynomial in pair
COEF_TYPE coeff(DEG_TYPE k) const
the th coefficient
uint64_t * serialized(uint64_t &size)
Dense_Univariate_Integer_Polynomial * hilbert_numerator_bigatti(const list< Monomial > &T, const WT_TYPE *grading)
the Bigatti algorithm to compute the Hilbert numerator
quick-’n-dirty Dense_Univariate integer polynomial class
DEG_TYPE poly_sugar() const
the sugar itself
void reduce_over_basis(Mutable_Polynomial **sp, const T &G, int comm_id=0)
Reduce the polynomial r over the basis G.
virtual void at_generation_tasks()
hook called while first generating polynomial
Definition: strategies.hpp:88
StrategyFlags
flag indicating which strategy to use for computation
Definition: strategies.hpp:34
list< Abstract_Polynomial * > reduce_basis(list< Abstract_Polynomial *>G)
Remove redundant polynomials from G.
SPolyCreationFlags
flag indicating which structure to use for an s-polynomial
DEG_TYPE size
number of slots for coefficients
Critical_Pair_XPlor(int i, Abstract_Polynomial *g, StrategyFlags strategy, vector< Abstract_Polynomial *>G)
create critical pair (f,g) where f is at index i
int qi
second polynomial in the critical pair
bool lcm_alike(const Monomial &t, const Monomial &u, const Critical_Pair_Basic *p)
Checks if the lcm of t and u is like the lcm stored in p.
NVAR_TYPE num_vars() const
number of variables
Definition: monomial.hpp:130
Polynomials that need arithmetic typically descend from this class.
Definition: polynomial.hpp:305
const Abstract_Polynomial * first() const
first polynomial in pair
virtual Monomial & leading_monomial() const override
Implementation of monomials.
Definition: monomial.hpp:69
interface to a monomial ordering
polynomial-related data for a sugar strategy
Encapsulates information about a polynomial ring for easy access: ground field, number of indetermina...
int pi
first polynomial in the critical pair
Critical_Pair_XPlor(int i, StrategyFlags strategy, Abstract_Polynomial *f)
create critical pair (f,0) where f is at index i
bool is_coprime(const Monomial &other) const
true iff this has no common factor with other.
Definition: monomial.cpp:243
DEG_TYPE degree(NVAR_TYPE i) const
Degree of th variable.
Definition: monomial.cpp:183
Pair_Strategy_Data * key
strategy used to sort critical pairs
void top_reduce(Mutable_Polynomial *s, Abstract_Polynomial *g, int comm_id)
reduce the polynomial **sp by *g
void set_strategy(Poly_Strategy_Data *psd)
sets the polynomial’s strategy to psd
Definition: polynomial.cpp:36
DEG_TYPE sugar()
returns sugar of this pair; use ONLY if with sugar strategy
void force_sugar(DEG_TYPE new_sugar)
for those times when a different sugar is appropriate
DEG_TYPE degree() const
the polynomial’s degree (exponent of largest nonzero term)
bool no_triplet(const T *p, const list< T *>C)
Checks whether p is in danger of forming a Buchberger triple with some pair listed in C...
quick-’n-dirty Dense_Univariate rational polynomial class
int proc
processor to which this pair has been assigned
MPZCOEF_TYPE denominator(DEG_TYPE k) const
returns the th denominator
DEG_TYPE degree() const
returns the polynomial’s degree
Controls the creation of s-polynomials.
unsigned ideal_dimension(NVAR_TYPE n, const Dense_Univariate_Integer_Polynomial *h1, const Dense_Univariate_Integer_Polynomial *h2)
computes the dimension of the ideal by subtracting the Hilbert numerators
Abstract_Polynomial * p
first polynomial in the critical pair
Abstract_Polynomial * find_reducer(Abstract_Polynomial *r, const T &G)
Find a polynomial in the basis G that can reduce r.
contains information on critical pairs by their index in the basis, in addition to the usual informat...
virtual Poly_Strategy_Data * strategy() const
strategy related information
Definition: polynomial.hpp:144