Gröbner basis project
Codebase for research into Gröbner basis computation
algorithm_buchberger_explorer_atlas.cpp
1 #ifndef __ALGORITHM_BUCHBERGER_EXPLORER_CPP_
2 #define __ALGORITHM_BUCHBERGER_EXPLORER_CPP_
3 
4 #include <vector>
5 using std::vector;
6 
7 #include "hilbert_functions.hpp"
8 
9 #include "algorithm_buchberger_basic.hpp"
10 #include "reduction_support.hpp"
11 
12 #include "algorithm_buchberger_explorer.hpp"
13 
14 #include "mpi.h"
15 
16 #define XPLOR_GENERATOR -1
17 #define PROC_UNASSIGNED -1
18 
25 public:
27 
29  Critical_Pair_XPlor(int i, unsigned strategy, Abstract_Polynomial * f)
30  : Critical_Pair_Basic(f, strategy)
31  {
32  pi = i; qi = XPLOR_GENERATOR; proc = PROC_UNASSIGNED;
33  }
36  int i, int j, unsigned strategy, vector<Abstract_Polynomial *>G
37  )
38  : Critical_Pair_Basic(G[i], G[j], strategy)
39  {
40  pi = i; qi = j; proc = PROC_UNASSIGNED;
41  }
44  int i, Abstract_Polynomial *g, unsigned strategy,
45  vector<Abstract_Polynomial *>G
46  )
47  : Critical_Pair_Basic(G[i], g, strategy)
48  {
49  pi = i; qi = G.size(); proc = PROC_UNASSIGNED;
50  }
52 
53 
55  int first_index() { return pi; }
57  int second_index() { return qi; }
59  DEG_TYPE sugar() {
60  return (static_cast<Pair_Sugar_Data *>(key))->pair_sugar();
61  }
63 
64 
66  void set_processor(int i) { proc = i; }
71  int get_processor() { return proc; }
73 protected:
75  int pi;
77  int qi;
79  int proc;
80 };
81 
88 bool second_HF_smaller(
93 ) {
94  // this adapts LessByHilbert in dynamic_engine.cpp
95  // in this case, we want the opposite result,
96  // so the return statement at the end is negated
97  bool result;
98  // first check the coefficients of the Hilbert polynomial
100  HPdiff -= *hp2;
101  if (not HPdiff.is_zero())
102  result = (HPdiff.numerator(HPdiff.degree()) < 0);
103  else // use Hilbert series
104  {
107  DEG_TYPE i = 0;
108  for ( /* already initialized */ ;
109  i <= h1->degree() and i <= h2->degree() and (*h1)[i] == (*h2)[i];
110  i++)
111  { /* taken care of in loop */ }
112  if (i > h1->degree())
113  {
114  if (i > h2->degree())
115  // the numerators are equal; second is not measurably better
116  result = false;
117  else
118  result = true;
119  }
120  else
121  {
122  if (i > h2->degree()) result = false;
123  else result = (*h1)[i] < (*h2)[i];
124  }
125  }
126  //cout << "\tfirst less than second? " << result << endl;
127  return not result;
128 }
129 
136 bool newer_HF_smaller(Monomial & t, unsigned i, Monomial & u, unsigned j,
139 ) {
140  // this adapts LessByHilbert in dynamic_engine.cpp
141  // in this case, we want the opposite result,
142  // so the return statement at the end is negated
143  bool result;
144  NVAR_TYPE n = t.num_vars();
145  // first check the coefficients of the Hilbert polynomial
146  Dense_Univariate_Rational_Polynomial HPdiff(*(HP[i]));
147  HPdiff -= *(HP[j]);
148  if (not HPdiff.is_zero())
149  result = (HPdiff.numerator(HPdiff.degree()) < 0);
150  else // use Hilbert series
151  {
154  DEG_TYPE i = 0;
155  for ( /* already initialized */ ;
156  i <= h1->degree() and i <= h2->degree() and (*h1)[i] == (*h2)[i];
157  i++)
158  { /* taken care of in loop */ }
159  if (i > h1->degree())
160  {
161  if (i > h2->degree())
162  { // the numerators are equal; break tie via lex
163  int i = 0;
164  while (i < n and t[i] == u[i]) ++i;
165  if (i == n) result = false;
166  else result = (t.degree(i) > u.degree(i));
167  }
168  else
169  result = true;
170  }
171  else
172  {
173  if (i > h2->degree()) result = false;
174  else result = (*h1)[i] < (*h2)[i];
175  }
176  }
177  //cout << "\tfirst less than second? " << result << endl;
178  return not result;
179 }
180 
194 void hilbert_chooser(
195  void * firstin, void * secondin, int * len, MPI_Datatype * dptr
196 ) {
197  int64_t * first = (int64_t *)firstin;
198  int64_t * second = (int64_t *)secondin;
199  bool undecided = true;
200  bool keep_second = true;
201  int64_t i;
202  int64_t n = first[1];
203  for (i = 2; undecided and i < n + 2; ++i)
204  undecided = keep_second = (first[i] < second[i]);
205  if (undecided)
206  for (/* i already assigned */ ; undecided and i < *len; ++i)
207  undecided = keep_second = (first[i] == second[i]);
208  // at this point, if we're still undecided then the functions are equal
209  if (not keep_second)
210  for (unsigned i = 0; i < *len; ++i)
211  second[i] = first[i];
212 }
213 
226  list<Critical_Pair_XPlor *> & P,
227  list<Critical_Pair_XPlor *> & Pass,
228  list<Critical_Pair_XPlor *> * Pcancel,
229  vector<Abstract_Polynomial *> & G,
231  unsigned strategy
232 ) {
233  //cout << "----------------------\n";
234  list<Critical_Pair_XPlor *> C;
235  //cout << "creating and pruning pairs, starting with:\n";
236  //cout << '\t' << P.size() << " pairs\n";
237  //cout << '\t' << G.size() << " polynomials\n";
238  // critical pairs with new polynomial
239  unsigned m = G.size();
240  for (unsigned i = 0; i < G.size(); ++i)
241  C.push_back(new Critical_Pair_XPlor(i, r, strategy, G));
242  //cout << "Created " << C.size() << " pairs with new polynomial\n";
243  // apply Buchberger's lcm criterion to new pairs
244  list<Critical_Pair_XPlor *> D;
245  while (C.size() != 0) {
246  Critical_Pair_XPlor * p = C.front();
247  C.pop_front();
248  if ((p->first()->leading_monomial().is_coprime(
249  p->second()->leading_monomial()))
250  or (no_triplet(p, C) and no_triplet(p, D))
251  )
252  D.push_back(p);
253  else {
254  //cout << "triplet prunes " << *p << endl;
255  delete p;
256  }
257  }
258  //cout << "After applying Buchberger's lcm criterion to new pairs, now have "
259  // << D.size() << " new pairs\n";
260  // apply Buchberger's gcd criterion
261  list<Critical_Pair_XPlor *> E;
262  while (D.size() != 0) {
263  Critical_Pair_XPlor * p = D.front();
264  D.pop_front();
265  if (!(p->first()->leading_monomial().is_coprime(
266  p->second()->leading_monomial())))
267  E.push_back(p);
268  else {
269  //cout << "gcd prunes " << *p << endl;
270  delete p;
271  }
272  }
273  //cout << "After applying Buchberger's gcd criterion to new pairs, now have "
274  // << E.size() << " new pairs\n";
275  // apply Buchberger's lcm criterion to old pairs
276  list<Critical_Pair_XPlor *> Q;
277  while (P.size() != 0) {
278  Critical_Pair_XPlor * p = P.front();
279  P.pop_front();
280  if (!(r->leading_monomial() | p->lcm())
283  )
284  Q.push_back(p);
285  else {
286  //cout << "triplet prunes " << *p << endl;
287  delete p;
288  }
289  }
290  //cout << "After applying Buchberger's lcm criterion to new pairs, now have "
291  // << Q.size() << " old pairs\n";
292  P = Q;
293  list <Critical_Pair_XPlor *> R;
294  while (Pass.size() != 0) {
295  Critical_Pair_XPlor * p = Pass.front();
296  Pass.pop_front();
297  if (!(r->leading_monomial() | p->lcm())
300  )
301  R.push_back(p);
302  else {
303  //cout << "\ttriplet prunes " << *p << endl;
304  Pcancel[p->get_processor()].push_back(p);
305  }
306  }
307  Pass = R;
308  // add new pairs to old pairs
309  for (Critical_Pair_XPlor * e : E)
310  P.push_back(e);
311  /*cout << "All pairs:\n";
312  for (list<Critical_Pair_XPlor *>::iterator pi = P.begin(); pi != P.end(); ++pi)
313  cout << '\t' << **pi << endl;
314  cout << "----------------------\n";*/
315 }
316 
318 
322 typedef struct {
324  int first;
326  int second;
328  DEG_TYPE sugar;
330 
331 list<Constant_Polynomial *> buchberger_explorer(
332  const vector<Abstract_Polynomial *> &F,
333  int method,
334  unsigned strategy,
335  WT_TYPE * strategy_weights,
336  const int comm_id,
337  const int comm_size
338 ) {
339  double r_bcast_time = 0;
340  unsigned number_of_spolys = 0;
341  vector<Abstract_Polynomial *> G; // basis
342  list<Critical_Pair_XPlor *> P; // critical pairs
343  list<Critical_Pair_XPlor *> Pass, * Pdel; // critical pair assignments
344  list<Critical_Pair_XPlor *> R; // reduced polynomials
345  list<Constant_Polynomial *> B; // end result
346  Polynomial_Ring & Rx = (*(F.begin()))->base_ring();
347  Monomial_Ordering * mord = (*(F.begin()))->monomial_ordering();
348  // set up MPI_Datatype
349  MPI_Datatype pair_type, pair_basetypes[2];
350  MPI_Aint pair_offsets[2], extent, lb;
351  int pair_block_counts[2];
352  pair_offsets[0] = 0; pair_basetypes[0] = MPI_INT; pair_block_counts[0] = 2;
353  MPI_Type_get_extent(MPI_INT, &lb, &extent);
354  pair_offsets[1] = 2 * extent; pair_basetypes[1] = MPI_UNSIGNED_LONG_LONG;
355  pair_block_counts[1] = 1;
356  MPI_Type_create_struct(2, pair_block_counts, pair_offsets, pair_basetypes, &pair_type);
357  MPI_Type_commit(&pair_type);
358  // set up MPI_Op
359  MPI_Op hilbert_comparison_op;
360  MPI_Op_create(hilbert_chooser, false, &hilbert_comparison_op);
361  // set up basis with generators
362  vector<Abstract_Polynomial *>::const_iterator Fi = F.begin();
363  NVAR_TYPE n = (*Fi)->number_of_variables();
364  for (unsigned i = 0; i < F.size(); ++i, ++Fi)
365  {
366  Abstract_Polynomial * fo = *Fi;
368  f->set_strategy(new Poly_Sugar_Data(f));
369  if (f->strategy() != nullptr) { f->strategy()->at_generation_tasks(); }
370  if (comm_id == 0) // only control needs to take care of critical pairs
371  P.push_back(new Critical_Pair_XPlor(i, strategy, f));
372  }
373  // main loop
374  bool verbose = false;
375  bool very_verbose = false;
376  unsigned min_todo = 0;
377  if (comm_id == 0) {
378  Pdel = new list<Critical_Pair_XPlor *>[comm_size];
380  = new Dense_Univariate_Rational_Polynomial *[comm_size];
382  = new Dense_Univariate_Integer_Polynomial *[comm_size];
384  = new Dense_Univariate_Integer_Polynomial *[comm_size];
385  min_todo = P.size();
386  }
387  MPI_Bcast(&min_todo, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
388  list<Monomial> T;
390  // create, reduce s-polynomials
392  while (min_todo != 0) {
393  /*for (int i = 0; i < comm_size; ++i) {
394  if (comm_id == i) {
395  cout << comm_id << "'s pairs\n";
396  for (Critical_Pair_XPlor * p : R)
397  cout << '\t' << comm_id << ' ' << *p << endl;
398  }
399  MPI_Barrier(MPI_COMM_WORLD);
400  }*/
401  if (comm_id == 0) {
402  if (P.size() > 0)
404  //cout << "estimate " << min_todo << " steps left\n";
405  // first select comm_id pairs for reduction, send to coprocessors, reduce
406  int i = 1;
408  for (/* already initialized */; !P.empty() and i < comm_size; ++i) {
409  Critical_Pair_XPlor * p = P.front();
410  p_new.first = p->first_index(); p_new.second = p->second_index();
411  p_new.sugar = p->sugar();
412  //report_front_pair(p, strategy);
413  P.pop_front();
414  p->set_processor(i);
415  Pass.push_back(p);
416  //cout << comm_id << " sending "
417  // << p_new.first << ',' << p_new.second << ':' << p_new.sugar
418  // << " to " << i << endl;
419  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
420  }
421  for (/* already initialized */; P.empty() and i < comm_size; ++i) {
422  p_new.first = p_new.second = XPLOR_GENERATOR; p_new.sugar = 0;
423  //cout << comm_id << " sending "
424  // << p_new.first << ',' << p_new.second << ':' << p_new.sugar
425  // << " to " << i << endl;
426  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
427  }
428  if (!P.empty()) {
429  Critical_Pair_XPlor * p = P.front();
430  P.pop_front();
431  //report_front_pair(p, strategy);
432  p_in.first = p->first_index(); p_in.second = p->second_index();
433  p_in.sugar = p->sugar();
434  p->set_processor(comm_id);
435  Pass.push_back(p);
436  //cout << comm_id << " sending "
437  // << p_in.first << ',' << p_in.second << ':' << p_in.sugar
438  // << " to " << comm_id << endl;
439  } else {
440  p_in.first = p_in.second = XPLOR_GENERATOR;
441  }
442  }
443  //MPI_Barrier(MPI_COMM_WORLD);
444  if (comm_id != 0)
445  MPI_Recv(&p_in, 1, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
446  if (p_in.first != XPLOR_GENERATOR) {
447  if (p_in.second == XPLOR_GENERATOR)
448  //R.push_back(new Critical_Pair_XPlor(p_in.first, p_in.second, p_in.sugar, F));
449  R.push_back(new Critical_Pair_XPlor(p_in.first, SUGAR_STRATEGY, F[p_in.first]));
450  else
451  R.push_back(new Critical_Pair_XPlor(p_in.first, p_in.second, p_in.sugar, G));
452  }
453  //MPI_Barrier(MPI_COMM_WORLD);
454  for (Critical_Pair_XPlor * p : R) {
455  // make s-poly
456  if ((s = p->s_polynomial()) == nullptr) {
457  s = p->s_polynomial(method, strategy);
458  ++number_of_spolys;
459  }
460  //cout << comm_id << " reducing s-poly "
461  // << p->first_index() << ',' << p->second_index() << endl;
462  //double start_time = MPI_Wtime();
463  if (!s->is_zero())
464  reduce_over_basis<vector<Abstract_Polynomial *>>(&s, G, comm_id);
465  //r_bcast_time += MPI_Wtime() - start_time;
466  p->set_spoly(s);
467  }
468  /*for (int i = 0; i < comm_size; ++i) {
469  if (comm_id == i)
470  cout << comm_id << " finished reduction with " << R.size() << " polynomials\n";
471  MPI_Barrier(MPI_COMM_WORLD);
472  }*/
473  // create, compare Hilbert functions
474  double start_time = MPI_Wtime();
475  unsigned winning_index = R.size() + 1;
476  Monomial * wt = nullptr;
477  Dense_Univariate_Rational_Polynomial * wHP = nullptr;
478  Dense_Univariate_Integer_Polynomial * wHS = nullptr;
480  Dense_Univariate_Integer_Polynomial * hsn = nullptr;
482  // loop through critical pairs to find optimal HF
483  list<Critical_Pair_XPlor *>::iterator Rbest = R.end();
484  //double start_time = MPI_Wtime();
485  for (
486  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
487  Ri != R.end();
488  /* advance manually */
489  ) {
490  s = (*Ri)->s_polynomial();
491  if (s->is_zero()) {
492  //cout << '\t' << comm_id << ' ' << (*Ri)->first_index() << ','
493  // << (*Ri)->second_index() << ": reduced to zero\n";
494  Critical_Pair_XPlor * p = *Ri;
495  if (Ri == R.begin()) {
496  R.pop_front();
497  Ri = R.begin();
498  } else {
499  list<Critical_Pair_XPlor *>::iterator Rdel = Ri;
500  ++Ri;
501  R.erase(Rdel);
502  }
503  delete p;
504  } else {
505  T.push_back(s->leading_monomial());
506  unsigned n = T.front().num_vars();
508  hsn = hilbert_second_numerator(n, hn);
509  unsigned d = ideal_dimension(n, hn, hsn);
510  hp = hilbert_polynomial(n, d, T, hn, hsn);
511  T.pop_back();
512  if (wt == nullptr) {
513  wt = &(s->leading_monomial()); wHP = hp; wHS = hn; Rbest = Ri;
514  } else {
515  if (second_HF_smaller(wHP, wHS, hp, hn)) {
516  wt = &(s->leading_monomial()); wHP = hp; wHS = hn; Rbest = Ri;
517  }
518  }
519  ++Ri;
520  }
521  }
522  //r_bcast_time = MPI_Wtime() - start_time;
523  //MPI_Barrier(MPI_COMM_WORLD);
524  // move result of winning reduction to basis; return others to P
525  /*cout << comm_id << " about to send Hilbert data\n";
526  for (int i = 1; i < comm_size; ++i) {
527  if (comm_id == i and wt != nullptr) {
528  cout << comm_id << "'s Hilbert data for " << *wt << ":\n";
529  cout << '\t' << *wHP << endl;
530  cout << '\t' << *wHS << endl;
531  }
532  MPI_Barrier(MPI_COMM_WORLD);
533  }*/
534  //MPI_Barrier(MPI_COMM_WORLD);
535  int64_t size_data = -1;
536  if (comm_id != 0) {
537  if (wHP != nullptr) {
538  size_data = wHP->degree();
539  MPI_Send(&size_data, 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
540  MPI_Send(wHP->numerator_array(), size_data + 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
541  MPI_Send(wHP->denominator_array(), size_data + 1, MPI_UINT64_T, 0, 0, MPI_COMM_WORLD);
542  delete wHP;
543  size_data = wHS->degree();
544  MPI_Send(&size_data, 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
545  MPI_Send(wHS->coefficient_array(), size_data + 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
546  delete wHS;
547  } else {
548  size_data = -1;
549  MPI_Send(&size_data, 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
550  }
551  }
552  int winners_id = -1;
553  //MPI_Barrier(MPI_COMM_WORLD);
554  if (comm_id == 0) {
555  /*cout << "sorting results\n";
556  cout << "My Hilbert data:\n";
557  if (wt == nullptr)
558  cout << "\tnothing\n";
559  else {
560  cout << '\t' << *wt << endl;
561  cout << '\t' << *wHP << endl;
562  cout << '\t' << *wHS << endl;
563  }*/
564  if (wHP != nullptr) winners_id = 0;
565  for (int i = 1; i < comm_size; ++i) {
566  MPI_Recv(&size_data, 1, MPI_INT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
567  if (size_data >= 0) {
568  int64_t * nums = new int64_t[size_data + 1];
569  uint64_t * denoms = new uint64_t[size_data + 1];
570  MPI_Recv(nums, size_data + 1, MPI_INT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
571  MPI_Recv(denoms, size_data + 1, MPI_UINT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
572  /*cout << "Received from " << i << ":\n\t";
573  for (unsigned k = 0; k < size_data + 1; ++k)
574  cout << nums[k] << '/' << denoms[k] << " , ";
575  cout << endl;*/
577  = new Dense_Univariate_Rational_Polynomial(size_data, nums, denoms);
578  //cout << "Received from " << i << ": " << *hp_in << endl;
579  delete [] nums; delete [] denoms;
580  MPI_Recv(&size_data, 1, MPI_UINT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
581  int64_t * coefs = new int64_t[size_data + 1];
582  MPI_Recv(coefs, size_data + 1, MPI_INT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
584  = new Dense_Univariate_Integer_Polynomial(size_data, coefs);
585  //cout << "Received from " << i << ": " << *hn_in << endl;
586  delete [] coefs;
587  if (wHP == nullptr) { // process 0 has exhausted its pairs
588  winners_id = i; wHP = hp_in; wHS = hn_in;
589  }
590  else if (second_HF_smaller(wHP, wHS, hp_in, hn_in)) {
591  if (wHP != nullptr) delete wHP;
592  if (wHS != nullptr) delete wHS;
593  winners_id = i; wHP = hp_in; wHS = hn_in;
594  } else {
595  delete hp_in; delete hn_in;
596  }
597  }
598  /*cout << "Received from " << i << ":\n";
599  cout << '\t' << *hp_in << endl;
600  cout << '\t' << *hn_in << endl;*/
601  }
602  if (wHP == nullptr) winners_id = -1;
603  /*cout << winners_id << " has the best HF\n";
604  if (winners_id != -1) {
605  cout << '\t' << *wHP << endl;
606  cout << '\t' << *wHS << endl;
607  }*/
608  }
609  MPI_Bcast(&winners_id, 1, MPI_INT, 0, MPI_COMM_WORLD);
610  if (winners_id >= 0) {
611  uint64_t * r_bcast;
612  uint64_t r_bcast_size;
613  uint64_t r_bcast_sugar;
615  if (comm_id == winners_id) {
616  // broadcast winning polynomial
617  s = (*Rbest)->s_polynomial();
618  r = new Constant_Polynomial(*s);
619  //cout << "selected " << (*Rbest)->first_index()
620  // << ',' << (*Rbest)->second_index() << ": " << r->leading_monomial()
621  // << endl;
622  Poly_Sugar_Data * sd = static_cast<Poly_Sugar_Data *>(s->strategy());
623  r->set_strategy(sd);
624  r_bcast_sugar = sd->poly_sugar();
625  s->set_strategy(nullptr);
626  delete s;
627  Critical_Pair_XPlor * p = *Rbest;
628  R.erase(Rbest);
629  delete p;
630  r_bcast = r->serialized(r_bcast_size);
631  }
632  //MPI_Barrier(MPI_COMM_WORLD);
633  MPI_Bcast(&r_bcast_size, 1, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
634  if (comm_id != winners_id)
635  r_bcast = new uint64_t [r_bcast_size];
636  //double start_time = MPI_Wtime();
637  MPI_Bcast(r_bcast, r_bcast_size, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
638  //r_bcast_time += MPI_Wtime() - start_time;
639  MPI_Bcast(&r_bcast_sugar, 1, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
640  // other processors need to create a copy of the new polynomial
641  if (comm_id != winners_id) {
642  r_bcast_size /= n + 1; // adjust size from # of words to # of terms
643  r = new Constant_Polynomial(Rx, mord, r_bcast_size, r_bcast);
644  Poly_Sugar_Data * rd = new Poly_Sugar_Data(r);
645  rd->force_sugar(r_bcast_sugar);
646  r->set_strategy(rd);
647  }
648  r_bcast_time += MPI_Wtime() - start_time;
649  delete [] r_bcast;
650  if (comm_id == 0) {
651  //cout << "added " << G.size() << ": " << r->leading_monomial() << endl;
652  very_verbose = false;
653  if (very_verbose) { cout << "\tin full "; r->println(); }
654  very_verbose = false;
655  gm_update_explorer(P, Pass, Pdel, G, r, strategy);
656  // remove useless pairs that were sent out to processors
657  for (int i = 1; i < comm_size; ++i) {
658  int outgoing = Pdel[i].size();
659  MPI_Send(&outgoing, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
660  while (Pdel[i].size() > 0) {
661  Critical_Pair_XPlor * p = Pdel[i].front();
663  Pdel[i].pop_front();
664  p_new.first = p->first_index();
665  p_new.second = p->second_index();
666  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
667  delete p;
668  }
669  }
670  // now delete my own redundant pairs
671  while (Pdel[0].size() > 0) {
672  Critical_Pair_XPlor * p = Pdel[0].front();
673  Pdel[0].pop_front();
674  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
675  while (
676  Ri != R.end() and
677  ((*Ri)->first_index() != p->first_index() or
678  (*Ri)->second_index() != p->second_index())
679  )
680  ++Ri;
681  // if the polynomial reduced to zero, it will not be in R
682  if (Ri != R.end())
683  R.erase(Ri);
684  delete p;
685  }
686  } else {
687  int incoming;
688  MPI_Recv(&incoming, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
689  for (int j = 0; j < incoming; ++j) {
690  MPI_Recv(&p_in, 1, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
691  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
692  while (
693  Ri != R.end() and
694  ((*Ri)->first_index() != p_in.first or
695  (*Ri)->second_index() != p_in.second)
696  )
697  ++Ri;
698  // if Ri has reduced to zero, it will not be in R
699  if (Ri != R.end()) {
700  Critical_Pair_XPlor * p = *Ri;
701  R.erase(Ri);
702  delete p;
703  }
704  }
705  }
706  //MPI_Barrier(MPI_COMM_WORLD);
707  G.push_back(r);
708  T.push_back(r->leading_monomial());
709  }
710  /*for (int i = 0; i < comm_size; ++i) {
711  if (comm_id == i) {
712  cout << comm_id << "'s basis\n";
713  for (Abstract_Polynomial * g : G)
714  cout << '\t' << comm_id << ' ' << *g << endl;
715  }
716  MPI_Barrier(MPI_COMM_WORLD);
717  }*/
718  if (comm_id == 0) min_todo = P.size();
719  unsigned my_todo = R.size();
720  unsigned all_todo;
721  MPI_Reduce(&my_todo, &all_todo, 1, MPI_UNSIGNED, MPI_MAX, 0, MPI_COMM_WORLD);
722  if (comm_id == 0)
723  min_todo = min_todo < all_todo ? all_todo : min_todo;
724  MPI_Bcast(&min_todo, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
725  //cout << comm_id << " understands there to be " << min_todo << " pairs remaining\n";
726  /*for (int i = 0; i < comm_size; ++i) {
727  if (comm_id == i)
728  for (Critical_Pair_XPlor * p : R)
729  cout << comm_id << "has " << p->first_index() << ',' << p->second_index() << endl;
730  }*/
731  //MPI_Barrier(MPI_COMM_WORLD);
732  }
733  if (comm_id != 0) {
734  for (Abstract_Polynomial * g : G)
735  delete g;
736  }
737  MPI_Type_free(&pair_type);
738  cout << comm_id << " reduced " << number_of_spolys << endl;
739  MPI_Barrier(MPI_COMM_WORLD);
740  int total_spolys;
741  MPI_Reduce(&number_of_spolys, &total_spolys, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
742  double max_bcast_time;
743  MPI_Reduce(&r_bcast_time, &max_bcast_time, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
744  if (comm_id == 0) {
745  cout << total_spolys << " s-polynomials computed and reduced\n";
746  cout << max_bcast_time << " seconds spent in communication\n";
747  // cleanup
748  cout << G.size() << " polynomials before interreduction\n";
749  //check_correctness(G, strategy);
750  list<Abstract_Polynomial *> G_final;
751  for (Abstract_Polynomial * g : G)
752  G_final.push_back(g);
753  G_final = reduce_basis(G_final);
754  cout << G_final.size() << " polynomials after interreduction\n";
755  //set<Constant_Polynomial *, smaller_lm> B;
756  for (Abstract_Polynomial * g : G_final) {
757  B.push_back(new Constant_Polynomial(*g));
758  //if (F.find(g) == F.end()) delete g;
759  }
760  delete [] Pdel;
761  }
762  return B;
763 }
764 
765 #endif
Critical_Pair_XPlor(int i, Abstract_Polynomial *g, unsigned strategy, vector< Abstract_Polynomial *>G)
create critical pair (f,g) where f is at index i
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
const COEF_TYPE * coefficient_array()
all the coefficients
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, unsigned strategy, Abstract_Polynomial *f)
create critical pair (f,0) where f is at index i
const Abstract_Polynomial * second() const
second polynomial in pair
Mutable_Polynomial * s
S-polynomial.
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, unsigned strategy)
Implementation of Gebauer-Moeller algorithm, with XPLOR critical pairs. Based on description in Becke...
Critical_Pair_XPlor(int i, int j, unsigned strategy, vector< Abstract_Polynomial *>G)
create critical pair (f,g) where f, g are at indices i, j
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_processor(int i)
record that this pair is assigned to processor i
int second_index()
returns index of second polynomial in pair
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
const COEF_TYPE * numerator_array()
all the numerators
list< Abstract_Polynomial * > reduce_basis(list< Abstract_Polynomial *>G)
Remove redundant polynomials from G.
DEG_TYPE size
number of slots for coefficients
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 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
DEG_TYPE degree() const
returns the polynomial’s degree
Controls the creation of s-polynomials.
const UCOEF_TYPE * denominator_array()
all the denominators
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
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