Gröbner basis project
Codebase for research into Gröbner basis computation
algorithm_buchberger_explorer_serial.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 "hilbert_functions.hpp"
22 
23 #include "algorithm_buchberger_basic.hpp"
24 
25 #include "algorithm_buchberger_explorer.hpp"
26 
27 #include "reduction_support.hpp"
28 
35 bool newer_HF_smaller(Monomial & t, unsigned i, Monomial & u, unsigned j,
38 ) {
39  // this adapts LessByHilbert in dynamic_engine.cpp
40  // in this case, we want the opposite result,
41  // so the return statement at the end is negated
42  bool result;
43  NVAR_TYPE n = t.num_vars();
44  // first check the coefficients of the Hilbert polynomial
45  Dense_Univariate_Rational_Polynomial HPdiff(*(HP[i]));
46  HPdiff -= *(HP[j]);
47  if (not HPdiff.is_zero())
48  result = (HPdiff.numerator(HPdiff.degree()) < 0);
49  else // use Hilbert series
50  {
53  DEG_TYPE i = 0;
54  for ( /* already initialized */ ;
55  i <= h1->degree() and i <= h2->degree() and (*h1)[i] == (*h2)[i];
56  i++)
57  { /* taken care of in loop */ }
58  if (i > h1->degree())
59  {
60  if (i > h2->degree())
61  { // the numerators are equal; break tie via lex
62  int i = 0;
63  while (i < n and t[i] == u[i]) ++i;
64  if (i == n) result = false;
65  else result = (t.degree(i) > u.degree(i));
66  }
67  else
68  result = true;
69  }
70  else
71  {
72  if (i > h2->degree()) result = false;
73  else result = (*h1)[i] < (*h2)[i];
74  }
75  }
76  //cout << "\tfirst less than second? " << result << endl;
77  return not result;
78 }
79 
81 
82 list<Constant_Polynomial *> buchberger_explorer(
83  const list<Abstract_Polynomial *> &F,
84  SPolyCreationFlags method,
85  StrategyFlags strategy,
86  WT_TYPE * strategy_weights,
87  const unsigned number_to_advance
88 ) {
89  unsigned number_of_spolys = 0;
90  list<Abstract_Polynomial *> G;
91  list<Critical_Pair_Basic *> P;
92  // set up basis with generators
93  for (Abstract_Polynomial * fo : F)
94  {
96  switch(strategy) {
97  case StrategyFlags::NORMAL_STRATEGY: break; // don't need polynomial data
98  case StrategyFlags::SUGAR_STRATEGY:
99  f->set_strategy(new Poly_Sugar_Data(f));
100  break;
101  case StrategyFlags::WSUGAR_STRATEGY:
102  f->set_strategy(new Poly_WSugar_Data(f, strategy_weights));
103  break;
104  default: break; // assume normal strategy
105  }
106  if (f->strategy() != nullptr) { f->strategy()->at_generation_tasks(); }
107  P.push_back(new Critical_Pair_Basic(f, strategy));
108  }
109  // main loop
110  bool verbose = false;
111  bool very_verbose = false;
112  Critical_Pair_Basic ** Pp = new Critical_Pair_Basic_Ptr[number_to_advance];
114  = new Dense_Univariate_Rational_Polynomial *[number_to_advance];
116  = new Dense_Univariate_Integer_Polynomial *[number_to_advance];
118  = new Dense_Univariate_Integer_Polynomial *[number_to_advance];
119  list<Monomial> T;
120  Mutable_Polynomial * s;
121  // create, reduce s-polynomials
122  while (!P.empty()) {
125  // first select number_to_advance pairs for reduction, store in Pp, reduce
126  unsigned i = 0;
127  for (/* already initialized */; !P.empty() and i < number_to_advance; ++i) {
128  Pp[i] = P.front();
129  report_front_pair(Pp[i], strategy);
130  P.pop_front();
131  // make s-poly
132  if ((s = Pp[i]->s_polynomial()) == nullptr)
133  s = Pp[i]->s_polynomial(method, strategy);
134  // cout << "Reducing s-poly "; s->println();
135  if (!s->is_zero())
136  reduce_over_basis(&s, G);
137  Pp[i]->set_spoly(s);
138  }
139  // create, compare Hilbert functions
140  while (i < number_to_advance) { Pp[i] = nullptr; ++i; }
141  unsigned winning_index = number_to_advance + 1;
142  for (i = 0; i < number_to_advance and Pp[i] != nullptr; ++i) {
143  if (Pp[i] != nullptr) {
144  s = Pp[i]->s_polynomial();
145  if (s->is_zero()) {
146  cout << "\treduced to zero\n";
147  HP[i] = nullptr;
148  HFN[i] = HSN[i]= nullptr;
149  } else {
150  T.push_back(s->leading_monomial());
151  unsigned n = T.front().num_vars();
152  HFN[i] = hilbert_numerator_bigatti(T);
153  HSN[i] = hilbert_second_numerator(n, HFN[i]);
154  unsigned d = ideal_dimension(n, HFN[i], HSN[i]);
155  HP[i] = hilbert_polynomial(n, d, T, HFN[i], HSN[i]);
156  T.pop_back();
157  if (winning_index > number_to_advance)
158  winning_index = i;
159  else {
160  if (newer_HF_smaller(
161  Pp[winning_index]->s_polynomial()->leading_monomial(),
162  winning_index,
163  Pp[i]->s_polynomial()->leading_monomial(),
164  i,
165  HP, HFN))
166  winning_index = i;
167  }
168  }
169  }
170  }
171  // move result of winning reduction to basis; return others to P
172  cout << "sorting results\n";
173  for (i = 0; i < number_to_advance and Pp[i] != nullptr; ++i) {
174  delete HP[i]; delete HFN[i]; delete HSN[i];
175  if (i != winning_index or Pp[i]->s_polynomial()->is_zero()) {
176  // not the winning pair: move to front of list of pairs
177  // (for further reduction)
178  if (not Pp[i]->s_polynomial()->is_zero())
179  P.push_front(Pp[i]);
180  else {
181  ++number_of_spolys;
182  cout << "\treduced to zero\n";
183  s = Pp[i]->s_polynomial();
184  //cout << "deleting " << s << endl;
185  delete s;
186  delete Pp[i];
187  }
188  } else { // winning pair
189  ++number_of_spolys;
190  cout << "selected "
191  << Pp[i]->first()->leading_monomial() << ','
192  << ((Pp[i]->is_generator()) ? 0 : Pp[i]->second()->leading_monomial())
193  << endl;
194  if (Pp[i]->is_generator()) delete Pp[i]->first();
195  s = Pp[i]->s_polynomial();
197  delete Pp[i];
198  // move strategy from s to r
199  r->set_strategy(s->strategy());
200  s->set_strategy(nullptr);
201  //cout << "deleting " << s << endl;
202  delete s;
203  cout << "\tadded " << r->leading_monomial() << endl;
204  if (very_verbose) { cout << "\tadded "; r->println(); }
205  gm_update(P, G, r, strategy);
206  }
207  }
208  }
209  delete [] Pp;
210  delete [] HP; delete [] HFN; delete [] HSN;
211  cout << number_of_spolys << " s-polynomials computed and reduced\n";
212  // cleanup
213  cout << G.size() << " polynomials before interreduction\n";
214  //check_correctness(G, strategy);
215  G = reduce_basis(G);
216  cout << G.size() << " polynomials after interreduction\n";
217  //set<Constant_Polynomial *, smaller_lm> B;
218  list<Constant_Polynomial *> B;
219  for (Abstract_Polynomial * g : G) {
220  B.push_back(new Constant_Polynomial(*g));
221  //if (F.find(g) == F.end()) delete g;
222  }
223  return B;
224 }
225 
226 #endif
The general class of a polynomial.
Definition: polynomial.hpp:101
virtual bool is_zero() const =0
is this polynomial zero?
A Constant_Polynomial is a polynomial that should not change.
virtual void set_spoly(Mutable_Polynomial *p)
sets the s-polynomial to p, for explorer methods
void sort_pairs_by_strategy(list< T *> &P)
Applies the strategy to find the “smallest” critical 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
void gm_update(list< Critical_Pair_Basic *> &P, list< Abstract_Polynomial *> &G, Abstract_Polynomial *r, StrategyFlags strategy)
Implementation of Gebauer-Moeller algorithm. Based on description in Becker and Weispfenning (1993)...
list< Constant_Polynomial * > buchberger_explorer(const list< Abstract_Polynomial *> &F, SPolyCreationFlags method, StrategyFlags strategy, WT_TYPE *strategy_weights, const unsigned number_to_advance)
Alternate implementation of Buchberger’s algorithm, for parallelization.
virtual Monomial & leading_monomial() const =0
leading monomial – call after sort_by_order()!
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 )
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
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
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
polynomial-related data for a weighted sugar strategy
Implementation of monomials.
Definition: monomial.hpp:69
polynomial-related data for a sugar strategy
virtual Mutable_Polynomial * s_polynomial()
to use only if s-polynomial is already computed by another method
DEG_TYPE degree(NVAR_TYPE i) const
Degree of th variable.
Definition: monomial.cpp:183
void set_strategy(Poly_Strategy_Data *psd)
sets the polynomial’s strategy to psd
Definition: polynomial.cpp:36
DEG_TYPE degree() const
the polynomial’s degree (exponent of largest nonzero term)
quick-’n-dirty Dense_Univariate rational polynomial class
Controls the creation of s-polynomials.
void report_critical_pairs(const list< T *>P, bool verbose=false)
A brief report on the number of critical pairs. If verbose is true, also lists them.
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
bool is_generator() const
whether this pair is from a generator
virtual Poly_Strategy_Data * strategy() const
strategy related information
Definition: polynomial.hpp:144