Gröbner basis project
Codebase for research into Gröbner basis computation
algorithm_buchberger_explorer_precustomreduce.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 
192 void gm_update_explorer(
193  list<Critical_Pair_XPlor *> & P,
194  list<Critical_Pair_XPlor *> & Pass,
195  list<Critical_Pair_XPlor *> * Pcancel,
196  vector<Abstract_Polynomial *> & G,
198  unsigned strategy
199 ) {
200  //cout << "----------------------\n";
201  list<Critical_Pair_XPlor *> C;
202  //cout << "creating and pruning pairs, starting with:\n";
203  //cout << '\t' << P.size() << " pairs\n";
204  //cout << '\t' << G.size() << " polynomials\n";
205  // critical pairs with new polynomial
206  unsigned m = G.size();
207  for (unsigned i = 0; i < G.size(); ++i)
208  C.push_back(new Critical_Pair_XPlor(i, r, strategy, G));
209  //cout << "Created " << C.size() << " pairs with new polynomial\n";
210  // apply Buchberger's lcm criterion to new pairs
211  list<Critical_Pair_XPlor *> D;
212  while (C.size() != 0) {
213  Critical_Pair_XPlor * p = C.front();
214  C.pop_front();
215  if ((p->first()->leading_monomial().is_coprime(
216  p->second()->leading_monomial()))
217  or (no_triplet(p, C) and no_triplet(p, D))
218  )
219  D.push_back(p);
220  else {
221  //cout << "triplet prunes " << *p << endl;
222  delete p;
223  }
224  }
225  //cout << "After applying Buchberger's lcm criterion to new pairs, now have "
226  // << D.size() << " new pairs\n";
227  // apply Buchberger's gcd criterion
228  list<Critical_Pair_XPlor *> E;
229  while (D.size() != 0) {
230  Critical_Pair_XPlor * p = D.front();
231  D.pop_front();
232  if (!(p->first()->leading_monomial().is_coprime(
233  p->second()->leading_monomial())))
234  E.push_back(p);
235  else {
236  //cout << "gcd prunes " << *p << endl;
237  delete p;
238  }
239  }
240  //cout << "After applying Buchberger's gcd criterion to new pairs, now have "
241  // << E.size() << " new pairs\n";
242  // apply Buchberger's lcm criterion to old pairs
243  list<Critical_Pair_XPlor *> Q;
244  while (P.size() != 0) {
245  Critical_Pair_XPlor * p = P.front();
246  P.pop_front();
247  if (!(r->leading_monomial() | p->lcm())
250  )
251  Q.push_back(p);
252  else {
253  //cout << "triplet prunes " << *p << endl;
254  delete p;
255  }
256  }
257  //cout << "After applying Buchberger's lcm criterion to new pairs, now have "
258  // << Q.size() << " old pairs\n";
259  P = Q;
260  list <Critical_Pair_XPlor *> R;
261  while (Pass.size() != 0) {
262  Critical_Pair_XPlor * p = Pass.front();
263  Pass.pop_front();
264  if (!(r->leading_monomial() | p->lcm())
267  )
268  R.push_back(p);
269  else {
270  //cout << "\ttriplet prunes " << *p << endl;
271  Pcancel[p->get_processor()].push_back(p);
272  }
273  }
274  Pass = R;
275  // add new pairs to old pairs
276  for (Critical_Pair_XPlor * e : E)
277  P.push_back(e);
278  /*cout << "All pairs:\n";
279  for (list<Critical_Pair_XPlor *>::iterator pi = P.begin(); pi != P.end(); ++pi)
280  cout << '\t' << **pi << endl;
281  cout << "----------------------\n";*/
282 }
283 
285 
289 typedef struct {
291  int first;
293  int second;
295  DEG_TYPE sugar;
297 
298 list<Constant_Polynomial *> buchberger_explorer(
299  const vector<Abstract_Polynomial *> &F,
300  int method,
301  unsigned strategy,
302  WT_TYPE * strategy_weights,
303  const int comm_id,
304  const int comm_size
305 ) {
306  double r_bcast_time = 0;
307  unsigned number_of_spolys = 0;
308  vector<Abstract_Polynomial *> G; // basis
309  list<Critical_Pair_XPlor *> P; // critical pairs
310  list<Critical_Pair_XPlor *> Pass, * Pdel; // critical pair assignments
311  list<Critical_Pair_XPlor *> R; // reduced polynomials
312  list<Constant_Polynomial *> B; // end result
313  Polynomial_Ring & Rx = (*(F.begin()))->base_ring();
314  Monomial_Ordering * mord = (*(F.begin()))->monomial_ordering();
315  // set up MPI_Datatype
316  MPI_Datatype pair_type, pair_basetypes[2];
317  MPI_Aint pair_offsets[2], extent, lb;
318  int pair_block_counts[2];
319  pair_offsets[0] = 0; pair_basetypes[0] = MPI_INT; pair_block_counts[0] = 2;
320  MPI_Type_get_extent(MPI_INT, &lb, &extent);
321  pair_offsets[1] = 2 * extent; pair_basetypes[1] = MPI_UNSIGNED_LONG_LONG;
322  pair_block_counts[1] = 1;
323  MPI_Type_create_struct(2, pair_block_counts, pair_offsets, pair_basetypes, &pair_type);
324  MPI_Type_commit(&pair_type);
325  // set up basis with generators
326  vector<Abstract_Polynomial *>::const_iterator Fi = F.begin();
327  NVAR_TYPE n = (*Fi)->number_of_variables();
328  for (unsigned i = 0; i < F.size(); ++i, ++Fi)
329  {
330  Abstract_Polynomial * fo = *Fi;
332  f->set_strategy(new Poly_Sugar_Data(f));
333  if (f->strategy() != nullptr) { f->strategy()->at_generation_tasks(); }
334  if (comm_id == 0) // only control needs to take care of critical pairs
335  P.push_back(new Critical_Pair_XPlor(i, strategy, f));
336  }
337  // main loop
338  bool verbose = false;
339  bool very_verbose = false;
340  unsigned min_todo = 0;
341  if (comm_id == 0) {
342  Pdel = new list<Critical_Pair_XPlor *>[comm_size];
344  = new Dense_Univariate_Rational_Polynomial *[comm_size];
346  = new Dense_Univariate_Integer_Polynomial *[comm_size];
348  = new Dense_Univariate_Integer_Polynomial *[comm_size];
349  min_todo = P.size();
350  }
351  MPI_Bcast(&min_todo, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
352  list<Monomial> T;
354  // create, reduce s-polynomials
356  while (min_todo != 0) {
357  /*for (int i = 0; i < comm_size; ++i) {
358  if (comm_id == i) {
359  cout << comm_id << "'s pairs\n";
360  for (Critical_Pair_XPlor * p : R)
361  cout << '\t' << comm_id << ' ' << *p << endl;
362  }
363  MPI_Barrier(MPI_COMM_WORLD);
364  }*/
365  //double start_time = MPI_Wtime();
366  if (comm_id == 0) {
367  if (P.size() > 0)
369  //cout << "estimate " << min_todo << " steps left\n";
370  // first select comm_id pairs for reduction, send to coprocessors, reduce
371  int i = 1;
373  for (/* already initialized */; !P.empty() and i < comm_size; ++i) {
374  Critical_Pair_XPlor * p = P.front();
375  p_new.first = p->first_index(); p_new.second = p->second_index();
376  p_new.sugar = p->sugar();
377  //report_front_pair(p, strategy);
378  P.pop_front();
379  p->set_processor(i);
380  Pass.push_back(p);
381  //cout << comm_id << " sending "
382  // << p_new.first << ',' << p_new.second << ':' << p_new.sugar
383  // << " to " << i << endl;
384  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
385  }
386  for (/* already initialized */; P.empty() and i < comm_size; ++i) {
387  p_new.first = p_new.second = XPLOR_GENERATOR; p_new.sugar = 0;
388  //cout << comm_id << " sending "
389  // << p_new.first << ',' << p_new.second << ':' << p_new.sugar
390  // << " to " << i << endl;
391  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
392  }
393  if (!P.empty()) {
394  Critical_Pair_XPlor * p = P.front();
395  P.pop_front();
396  //report_front_pair(p, strategy);
397  p_in.first = p->first_index(); p_in.second = p->second_index();
398  p_in.sugar = p->sugar();
399  p->set_processor(comm_id);
400  Pass.push_back(p);
401  //cout << comm_id << " sending "
402  // << p_in.first << ',' << p_in.second << ':' << p_in.sugar
403  // << " to " << comm_id << endl;
404  } else {
405  p_in.first = p_in.second = XPLOR_GENERATOR;
406  }
407  }
408  //MPI_Barrier(MPI_COMM_WORLD);
409  if (comm_id != 0)
410  MPI_Recv(&p_in, 1, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
411  if (p_in.first != XPLOR_GENERATOR) {
412  if (p_in.second == XPLOR_GENERATOR)
413  //R.push_back(new Critical_Pair_XPlor(p_in.first, p_in.second, p_in.sugar, F));
414  R.push_back(new Critical_Pair_XPlor(p_in.first, SUGAR_STRATEGY, F[p_in.first]));
415  else
416  R.push_back(new Critical_Pair_XPlor(p_in.first, p_in.second, SUGAR_STRATEGY, G));
417  }
418  //r_bcast_time += MPI_Wtime() - start_time;
419  //MPI_Barrier(MPI_COMM_WORLD);
420  double start_time = MPI_Wtime();
421  for (Critical_Pair_XPlor * p : R) {
422  // make s-poly
423  if ((s = p->s_polynomial()) == nullptr) {
424  s = p->s_polynomial(method, strategy);
425  ++number_of_spolys;
426  }
427  //cout << comm_id << " reducing s-poly "
428  // << p->first_index() << ',' << p->second_index() << endl;
429  //double start_time = MPI_Wtime();
430  if (!s->is_zero())
431  reduce_over_basis<vector<Abstract_Polynomial *>>(&s, G, comm_id);
432  //r_bcast_time += MPI_Wtime() - start_time;
433  p->set_spoly(s);
434  }
435  r_bcast_time += MPI_Wtime() - start_time;
436  /*for (int i = 0; i < comm_size; ++i) {
437  if (comm_id == i)
438  cout << comm_id << " finished reduction with " << R.size() << " polynomials\n";
439  MPI_Barrier(MPI_COMM_WORLD);
440  }*/
441  // create, compare Hilbert functions
442  //MPI_Barrier(MPI_COMM_WORLD);
443  //double start_time = MPI_Wtime();
444  unsigned winning_index = R.size() + 1;
445  Monomial * wt = nullptr;
446  Dense_Univariate_Rational_Polynomial * wHP = nullptr;
447  Dense_Univariate_Integer_Polynomial * wHS = nullptr;
449  Dense_Univariate_Integer_Polynomial * hsn = nullptr;
451  // loop through critical pairs to find optimal HF
452  list<Critical_Pair_XPlor *>::iterator Rbest = R.end();
453  //double start_time = MPI_Wtime();
454  for (
455  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
456  Ri != R.end();
457  /* advance manually */
458  ) {
459  s = (*Ri)->s_polynomial();
460  if (s->is_zero()) {
461  //cout << '\t' << comm_id << ' ' << (*Ri)->first_index() << ','
462  // << (*Ri)->second_index() << ": reduced to zero\n";
463  Critical_Pair_XPlor * p = *Ri;
464  if (Ri == R.begin()) {
465  R.pop_front();
466  Ri = R.begin();
467  } else {
468  list<Critical_Pair_XPlor *>::iterator Rdel = Ri;
469  ++Ri;
470  R.erase(Rdel);
471  }
472  delete p;
473  } else {
474  T.push_back(s->leading_monomial());
475  unsigned n = T.front().num_vars();
477  hsn = hilbert_second_numerator(n, hn);
478  unsigned d = ideal_dimension(n, hn, hsn);
479  hp = hilbert_polynomial(n, d, T, hn, hsn);
480  T.pop_back();
481  if (wt == nullptr) {
482  wt = &(s->leading_monomial()); wHP = hp; wHS = hn; Rbest = Ri;
483  } else {
484  if (second_HF_smaller(wHP, wHS, hp, hn)) {
485  wt = &(s->leading_monomial()); wHP = hp; wHS = hn; Rbest = Ri;
486  }
487  }
488  ++Ri;
489  }
490  }
491  //r_bcast_time = MPI_Wtime() - start_time;
492  //MPI_Barrier(MPI_COMM_WORLD);
493  // move result of winning reduction to basis; return others to P
494  /*cout << comm_id << " about to send Hilbert data\n";
495  for (int i = 1; i < comm_size; ++i) {
496  if (comm_id == i and wt != nullptr) {
497  cout << comm_id << "'s Hilbert data for " << *wt << ":\n";
498  cout << '\t' << *wHP << endl;
499  cout << '\t' << *wHS << endl;
500  }
501  MPI_Barrier(MPI_COMM_WORLD);
502  }*/
503  //MPI_Barrier(MPI_COMM_WORLD);
504  //double start_time = MPI_Wtime();
505  int64_t size_data = -1;
506  if (comm_id != 0) {
507  if (wHP != nullptr) {
508  size_data = wHP->degree();
509  MPI_Send(&size_data, 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
510  MPI_Send(wHP->numerator_array(), size_data + 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
511  MPI_Send(wHP->denominator_array(), size_data + 1, MPI_UINT64_T, 0, 0, MPI_COMM_WORLD);
512  delete wHP;
513  size_data = wHS->degree();
514  MPI_Send(&size_data, 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
515  MPI_Send(wHS->coefficient_array(), size_data + 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
516  delete wHS;
517  } else {
518  size_data = -1;
519  MPI_Send(&size_data, 1, MPI_INT64_T, 0, 0, MPI_COMM_WORLD);
520  }
521  }
522  int winners_id = -1;
523  //MPI_Barrier(MPI_COMM_WORLD);
524  if (comm_id == 0) {
525  /*cout << "sorting results\n";
526  cout << "My Hilbert data:\n";
527  if (wt == nullptr)
528  cout << "\tnothing\n";
529  else {
530  cout << '\t' << *wt << endl;
531  cout << '\t' << *wHP << endl;
532  cout << '\t' << *wHS << endl;
533  }*/
534  if (wHP != nullptr) winners_id = 0;
535  for (int i = 1; i < comm_size; ++i) {
536  int64_t size_data;
537  MPI_Recv(&size_data, 1, MPI_INT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
538  if (size_data >= 0) {
539  int64_t * nums = new int64_t[size_data + 1];
540  uint64_t * denoms = new uint64_t[size_data + 1];
541  MPI_Recv(nums, size_data + 1, MPI_INT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
542  MPI_Recv(denoms, size_data + 1, MPI_UINT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
543  /*cout << "Received from " << i << ":\n\t";
544  for (unsigned k = 0; k < size_data + 1; ++k)
545  cout << nums[k] << '/' << denoms[k] << " , ";
546  cout << endl;*/
548  = new Dense_Univariate_Rational_Polynomial(size_data, nums, denoms);
549  //cout << "Received from " << i << ": " << *hp_in << endl;
550  delete [] nums; delete [] denoms;
551  MPI_Recv(&size_data, 1, MPI_UINT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
552  int64_t * coefs = new int64_t[size_data + 1];
553  MPI_Recv(coefs, size_data + 1, MPI_INT64_T, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
555  = new Dense_Univariate_Integer_Polynomial(size_data, coefs);
556  //cout << "Received from " << i << ": " << *hn_in << endl;
557  delete [] coefs;
558  if (wHP == nullptr) { // process 0 has exhausted its pairs
559  winners_id = i; wHP = hp_in; wHS = hn_in;
560  }
561  else if (second_HF_smaller(wHP, wHS, hp_in, hn_in)) {
562  if (wHP != nullptr) delete wHP;
563  if (wHS != nullptr) delete wHS;
564  winners_id = i; wHP = hp_in; wHS = hn_in;
565  } else {
566  delete hp_in; delete hn_in;
567  }
568  }
569  /*cout << "Received from " << i << ":\n";
570  cout << '\t' << *hp_in << endl;
571  cout << '\t' << *hn_in << endl;*/
572  }
573  if (wHP == nullptr) winners_id = -1;
574  /*cout << winners_id << " has the best HF\n";
575  if (winners_id != -1) {
576  cout << '\t' << *wHP << endl;
577  cout << '\t' << *wHS << endl;
578  }*/
579  }
580  MPI_Bcast(&winners_id, 1, MPI_INT, 0, MPI_COMM_WORLD);
581  if (winners_id >= 0) {
582  uint64_t * r_bcast;
583  uint64_t r_bcast_size;
584  uint64_t r_bcast_sugar;
586  if (comm_id == winners_id) {
587  // broadcast winning polynomial
588  s = (*Rbest)->s_polynomial();
589  r = new Constant_Polynomial(*s);
590  //cout << "selected " << (*Rbest)->first_index()
591  // << ',' << (*Rbest)->second_index() << ": " << r->leading_monomial()
592  // << endl;
593  Poly_Sugar_Data * sd = static_cast<Poly_Sugar_Data *>(s->strategy());
594  r->set_strategy(sd);
595  r_bcast_sugar = sd->poly_sugar();
596  s->set_strategy(nullptr);
597  delete s;
598  Critical_Pair_XPlor * p = *Rbest;
599  R.erase(Rbest);
600  delete p;
601  r_bcast = r->serialized(r_bcast_size);
602  }
603  //MPI_Barrier(MPI_COMM_WORLD);
604  MPI_Bcast(&r_bcast_size, 1, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
605  if (comm_id != winners_id)
606  r_bcast = new uint64_t [r_bcast_size];
607  //double start_time = MPI_Wtime();
608  MPI_Bcast(r_bcast, r_bcast_size, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
609  //r_bcast_time += MPI_Wtime() - start_time;
610  MPI_Bcast(&r_bcast_sugar, 1, MPI_UINT64_T, winners_id, MPI_COMM_WORLD);
611  // other processors need to create a copy of the new polynomial
612  if (comm_id != winners_id) {
613  r_bcast_size /= n + 1; // adjust size from # of words to # of terms
614  r = new Constant_Polynomial(Rx, mord, r_bcast_size, r_bcast);
615  Poly_Sugar_Data * rd = new Poly_Sugar_Data(r);
616  rd->force_sugar(r_bcast_sugar);
617  r->set_strategy(rd);
618  }
619  //MPI_Barrier(MPI_COMM_WORLD);
620  //r_bcast_time += MPI_Wtime() - start_time;
621  delete [] r_bcast;
622  //double start_time = MPI_Wtime();
623  if (comm_id == 0) {
624  //cout << "added " << G.size() << ": " << r->leading_monomial() << endl;
625  very_verbose = false;
626  if (very_verbose) { cout << "\tin full "; r->println(); }
627  very_verbose = false;
628  gm_update_explorer(P, Pass, Pdel, G, r, strategy);
629  // remove useless pairs that were sent out to processors
630  for (int i = 1; i < comm_size; ++i) {
631  int outgoing = Pdel[i].size();
632  MPI_Send(&outgoing, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
633  while (Pdel[i].size() > 0) {
634  Critical_Pair_XPlor * p = Pdel[i].front();
636  Pdel[i].pop_front();
637  p_new.first = p->first_index();
638  p_new.second = p->second_index();
639  MPI_Send(&p_new, 1, pair_type, i, 0, MPI_COMM_WORLD);
640  delete p;
641  }
642  }
643  // now delete my own redundant pairs
644  while (Pdel[0].size() > 0) {
645  Critical_Pair_XPlor * p = Pdel[0].front();
646  Pdel[0].pop_front();
647  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
648  while (
649  Ri != R.end() and
650  ((*Ri)->first_index() != p->first_index() or
651  (*Ri)->second_index() != p->second_index())
652  )
653  ++Ri;
654  // if the polynomial reduced to zero, it will not be in R
655  if (Ri != R.end())
656  R.erase(Ri);
657  delete p;
658  }
659  } else {
660  int incoming;
661  MPI_Recv(&incoming, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
662  for (int j = 0; j < incoming; ++j) {
663  MPI_Recv(&p_in, 1, pair_type, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
664  list<Critical_Pair_XPlor *>::iterator Ri = R.begin();
665  while (
666  Ri != R.end() and
667  ((*Ri)->first_index() != p_in.first or
668  (*Ri)->second_index() != p_in.second)
669  )
670  ++Ri;
671  // if Ri has reduced to zero, it will not be in R
672  if (Ri != R.end()) {
673  Critical_Pair_XPlor * p = *Ri;
674  R.erase(Ri);
675  delete p;
676  }
677  }
678  }
679  //r_bcast_time += MPI_Wtime() - start_time;
680  //MPI_Barrier(MPI_COMM_WORLD);
681  G.push_back(r);
682  T.push_back(r->leading_monomial());
683  }
684  //r_bcast_time += MPI_Wtime() - start_time;
685  /*for (int i = 0; i < comm_size; ++i) {
686  if (comm_id == i) {
687  cout << comm_id << "'s basis\n";
688  for (Abstract_Polynomial * g : G)
689  cout << '\t' << comm_id << ' ' << *g << endl;
690  }
691  MPI_Barrier(MPI_COMM_WORLD);
692  }*/
693  //double start_time = MPI_Wtime();
694  if (comm_id == 0) min_todo = P.size();
695  unsigned my_todo = R.size();
696  unsigned all_todo;
697  MPI_Reduce(&my_todo, &all_todo, 1, MPI_UNSIGNED, MPI_MAX, 0, MPI_COMM_WORLD);
698  if (comm_id == 0)
699  min_todo = min_todo < all_todo ? all_todo : min_todo;
700  MPI_Bcast(&min_todo, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
701  //r_bcast_time += MPI_Wtime() - start_time;
702  //cout << comm_id << " understands there to be " << min_todo << " pairs remaining\n";
703  /*for (int i = 0; i < comm_size; ++i) {
704  if (comm_id == i)
705  for (Critical_Pair_XPlor * p : R)
706  cout << comm_id << "has " << p->first_index() << ',' << p->second_index() << endl;
707  MPI_Barrier(MPI_COMM_WORLD);
708  }*/
709  //MPI_Barrier(MPI_COMM_WORLD);
710  }
711  if (comm_id != 0) {
712  for (Abstract_Polynomial * g : G)
713  delete g;
714  }
715  MPI_Type_free(&pair_type);
716  cout << comm_id << " reduced " << number_of_spolys << endl;
717  MPI_Barrier(MPI_COMM_WORLD);
718  int total_spolys;
719  MPI_Reduce(&number_of_spolys, &total_spolys, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
720  double max_bcast_time;
721  MPI_Reduce(&r_bcast_time, &max_bcast_time, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
722  for (unsigned i = 0; i < comm_size; ++i) {
723  if (comm_id == i)
724  cout << comm_id << " took " << r_bcast_time << " milliseconds on timed segment(s).\n";
725  MPI_Barrier(MPI_COMM_WORLD);
726  }
727  if (comm_id == 0) {
728  cout << total_spolys << " s-polynomials computed and reduced\n";
729  cout << max_bcast_time << " seconds spent on timed segment(s)\n";
730  // cleanup
731  cout << G.size() << " polynomials before interreduction\n";
732  //check_correctness(G, strategy);
733  list<Abstract_Polynomial *> G_final;
734  for (Abstract_Polynomial * g : G)
735  G_final.push_back(g);
736  G_final = reduce_basis(G_final);
737  cout << G_final.size() << " polynomials after interreduction\n";
738  //set<Constant_Polynomial *, smaller_lm> B;
739  for (Abstract_Polynomial * g : G_final) {
740  B.push_back(new Constant_Polynomial(*g));
741  //if (F.find(g) == F.end()) delete g;
742  }
743  delete [] Pdel;
744  }
745  return B;
746 }
747 
748 #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
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.
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
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.
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