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