Gröbner basis project
Codebase for research into Gröbner basis computation
f4_reduction_serial.cpp
1 #ifndef __F4_REDUCTION_SERIAL_CPP__
2 #define __F4_REDUCTION_SERIAL_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 "f4_reduction_serial.hpp"
22 #include "algorithm_buchberger_basic.hpp"
23 
24 void find_position(
25  list<Monomial *>::iterator & ti,
26  list<Monomial *>::iterator stop,
27  list<Abstract_Polynomial *>::iterator & pi,
28  const Monomial & u
29 ) {
30  while (ti != stop and **ti > u) {
31  ++ti;
32  ++pi;
33  }
34 }
35 
37  const list<Critical_Pair_Basic *> & P,
38  const list<Abstract_Polynomial *> & B
39 ) : G(B), Rx(P.front()->first()->base_ring()) {
40  cols = 0;
41  A = nullptr;
42  R_build.clear();
43  M_build.clear();
44  strategies.clear();
45  auto mi = M_build.begin();
46  auto ri = R_build.begin();
47  mord = P.front()->first()->monomial_ordering();
48  for (auto p : P) {
49  strategies.push_back(new Poly_Sugar_Data(p->first()));
50  if (strategies.back() != nullptr) {
51  strategies.back()->at_generation_tasks(p->first_multiplier());
52  if (p->second() != nullptr)
53  strategies.back()->pre_reduction_tasks(p->second_multiplier().log(), *(p->second()));
54  }
55  add_monomials(mi, ri, p->first(), p->first_multiplier(), true);
56  if (p->second() != nullptr) {
57  *ri = const_cast<Abstract_Polynomial *>(p->second());
58  ++ri; ++mi;
59  add_monomials(mi, ri, p->second(), p->second_multiplier());
60  }
61  }
62  // for each monomial, find an appropriate reducer
63  mi = M_build.begin();
64  ri = R_build.begin();
65  while (mi != M_build.end()) {
66  auto g = G.begin();
67  bool found = ((*ri) != nullptr);
68  while (not found and g != G.end()) {
69  if ((**mi).divisible_by((*g)->leading_monomial())) {
70  found = true;
71  *ri = *g;
72  Monomial u(**mi);
73  u /= (*g)->leading_monomial();
74  add_monomials(mi, ri, *g, u);
75  g = G.end();
76  }
77  ++g;
78  }
79  ++mi;
80  ++ri;
81  }
82  initialize_many(P);
83 }
84 
85 void F4_Reduction_Data::initialize_many(const list<Critical_Pair_Basic *> & P) {
86  cols = M_build.size();
87  Abstract_Polynomial * p = const_cast<Abstract_Polynomial *>(P.front()->first());
88  Prime_Field & F = const_cast<Prime_Field &>(p->ground_field());
89  rows = P.size();
90  heads.resize(rows);
91  nonzero_entries.resize(rows);
92  for (unsigned i = 0; i < rows; ++i)
93  heads[i] = 0;
94  A = (Prime_Field_Element *)malloc(
95  sizeof(Prime_Field_Element)*cols*rows*P.size()
96  );
97  for (unsigned i = 0; i < cols*rows; ++i)
98  A[i].assign(0, &F);
99  M.clear();
100  for (auto mi = M_build.begin(); mi != M_build.end(); ++mi)
101  M.push_back(*mi);
102  strategies.resize(P.size());
103  for (unsigned i = 0; i < strategies.size(); ++i)
104  strategies[i] = nullptr;
105  unsigned row = 0;
106  for (auto cp : P) {
107  auto p = cp->first();
108  auto t = cp->first_multiplier();
109  auto pi = p->new_iterator();
110  nonzero_entries[row] = 0;
111  auto mi = M_build.begin();
112  bool head_found = false;
113  for (unsigned i = 0; i < cols; ++i) {
114  if ((not pi->fellOff()) and M[i]->like_multiple(pi->currMonomial(), t)) {
115  *(A + row*cols + i) = pi->currCoeff();
116  pi->moveRight();
117  ++nonzero_entries[row];
118  if (not head_found) {
119  heads[row] = i;
120  head_found = true;
121  }
122  }
123  ++mi;
124  }
125  delete pi;
126  strategies[row] = new Poly_Sugar_Data(cp->first());
127  delete cp;
128  ++row;
129  }
130  R.resize(R_build.size());
131  unsigned i = 0;
132  for (Abstract_Polynomial * r : R_build)
133  R[i++] = r;
134 }
135 
137  list<Monomial *>::iterator & t1,
138  list<Abstract_Polynomial *>::iterator & r1,
139  const Abstract_Polynomial *g,
140  const Monomial & u,
141  bool new_row
142 ) {
143  NVAR_TYPE n = g->leading_monomial().num_vars();
144  if (new_row) {
145  t1 = M_build.begin();
146  r1 = R_build.begin();
147  Monomial * t = new Monomial(g->leading_monomial());
148  (*t) *= u;
149  find_position(t1, M_build.end(), r1, *t);
150  if (t1 == M_build.end()) {
151  M_build.push_back(t);
152  R_build.push_back(nullptr);
153  t1 = M_build.begin();
154  r1 = R_build.begin();
155  auto t2(t1);
156  ++t2;
157  while (t2 != M_build.end()) {
158  ++t1; ++r1;
159  ++t2;
160  }
161  } else if (**t1 != *t) {
162  M_build.insert(t1, t);
163  R_build.insert(r1, nullptr);
164  --t1; --r1;
165  }
166  }
167  auto t2(t1); auto r2(r1);
168  Polynomial_Iterator * pi = g->new_iterator();
169  pi->moveRight();
170  while (not (pi->fellOff())) {
171  Monomial * t = new Monomial(pi->currMonomial());
172  (*t) *= u;
173  find_position(t2, M_build.end(), r2, *t);
174  if (t2 == M_build.end()) {
175  M_build.push_back(t);
176  R_build.push_back(nullptr);
177  } else if (**t2 != *t) {
178  M_build.insert(t2, t);
179  R_build.insert(r2, nullptr);
180  }
181  pi->moveRight();
182  }
183  delete pi;
184 }
185 
187  free(A);
188  for (auto strat : strategies) {
189  if (strat != nullptr)
190  delete strat;
191  }
192 }
193 
194 void F4_Reduction_Data::print_matrix(bool show_data) {
195  if (show_data) {
196  for (auto m : M)
197  cout << *m << ", ";
198  cout << endl;
199  }
200  for (unsigned i = 0; i < rows; ++i) {
201  cout << "( ";
202  for (unsigned j = 0; j < cols; ++j) {
203  cout << A[i*cols + j];
204  if (j < cols - 1)
205  cout << ", ";
206  }
207  cout << ")\n";
208  }
209 }
210 
212  for (unsigned i = 0; i < cols; ++i) {
213  cout << *(M[i]) << " to be reduced by ";
214  if (R[i] == nullptr)
215  cout << "none\n";
216  else
217  cout << R[i]->leading_monomial() << endl;
218  }
219 }
220 
222  bool is_zero_so_far = true;
223  for (unsigned i = 0; is_zero_so_far and i < rows; ++i)
224  is_zero_so_far = is_zero_so_far and (nonzero_entries[i] == 0);
225  return is_zero_so_far;
226 }
227 
229  NVAR_TYPE n = M[0]->num_vars();
230  EXP_TYPE * u = new EXP_TYPE[n]; // exponents of multiplier
231  // loop through rows
232  for (unsigned k = 0; k < rows; ++k) {
233  // loop through our terms
234  for (unsigned i = heads[k]; i < cols; ++i) {
235  // do we need to reduce this term, and can we?
236  if ((not A[k*cols + i].is_zero()) and R[i] != nullptr) {
237  // get reducer for this monomial
238  const Abstract_Polynomial * g = R[i];
239  Polynomial_Iterator * gi = g->new_iterator();
240  // determine multiplier
241  for (NVAR_TYPE k = 0; k < n; ++k)
242  u[k] = (*M[i])[k] - (gi->currMonomial())[k];
243  // determine reduction coefficient
244  Prime_Field_Element a(A[k*cols + i]*gi->currCoeff().inverse());
245  // loop through g's terms
246  // by construction, monomials of u*g[i] should already appear in M,
247  // so the line marked *** SHOULD NOT need a guard
248  unsigned j = i;
249  while (not gi->fellOff()) {
250  const Monomial & t = gi->currMonomial();
251  while (not M[j]->like_multiple(u, t)) ++j;
252  bool was_zero = A[k*cols + j].is_zero();
253  A[k*cols + j] -= a*gi->currCoeff();
254  if (was_zero and not A[k*cols + j].is_zero())
255  ++nonzero_entries[k];
256  else if (not was_zero and A[k*cols + j].is_zero())
257  --nonzero_entries[k];
258  ++j;
259  gi->moveRight();
260  }
261  if (strategies[k] != nullptr)
262  strategies[k]->pre_reduction_tasks(u, *g);
263  advance_head(k);
264  delete gi;
265  }
266  }
267  }
268  delete [] u;
269 }
270 
272  for (unsigned i = 0; i < rows; ++i) {
273  if (nonzero_entries[i] > 0) {
274  for (unsigned j = 0; j < rows; ++j) {
275  if (j != i) {
276  auto c = heads[i];
277  if (not A[j*cols + c].is_zero()) {
278  auto a = A[j*cols + c];
279  a *= -A[i*cols + c].inverse();
280  while (c < cols) {
281  if (not A[i*cols + c].is_zero()) {
282  bool was_zero = A[j*cols + c].is_zero();
283  A[j*cols + c] += a*A[i*cols + c];
284  if (was_zero and not A[j*cols + c].is_zero())
285  ++nonzero_entries[j];
286  else if (not was_zero and A[j*cols + c].is_zero())
287  --nonzero_entries[j];
288  }
289  ++c;
290  }
291  if (heads[i] == heads[j])
292  advance_head(j);
293  }
294  }
295  }
296  }
297  }
298 }
299 
300 vector<Constant_Polynomial *> F4_Reduction_Data::finalize() {
301  vector<Constant_Polynomial *> result;
302  NVAR_TYPE n = M[0]->num_vars();
303  for (unsigned i = 0; i < rows; ++i) {
304  if (nonzero_entries[i] == 0) {
305  delete strategies[i];
306  strategies[i] = nullptr;
307  } else {
308  Monomial * M_final = (Monomial *)malloc(sizeof(Monomial)*nonzero_entries[i]);
309  Prime_Field_Element * A_final
311  unsigned k = 0;
312  Prime_Field_Element scale(A[i*cols + heads[i]].inverse(), A[i*cols + heads[i]].field());
313  for (unsigned j = heads[i]; k < nonzero_entries[i] and j < cols; ++j) {
314  if (not A[i*cols + j].is_zero()) {
315  A_final[k] = A[i*cols + j]*scale;
316  M_final[k].common_initialization();
317  M_final[k].initialize_exponents(n);
318  M_final[k].set_monomial_ordering(mord);
319  M_final[k] = *M[j];
320  ++k;
321  }
322  }
323  result.push_back(new Constant_Polynomial(
324  nonzero_entries[i],
325  Rx,
326  M_final, A_final,
327  mord
328  ));
329  result.back()->set_strategy(strategies[i]);
330  strategies[i] = nullptr;
331  free(M_final);
332  free(A_final);
333  }
334  }
335  return result;
336 }
337 
338 list<Constant_Polynomial *> f4_control(const list<Abstract_Polynomial *> &F) {
339  unsigned number_of_spolys = 0;
340  list<Abstract_Polynomial *> G;
341  list<Critical_Pair_Basic *> P;
342  // set up basis with generators
343  for (Abstract_Polynomial * fo : F)
344  {
346  f->set_strategy(new Poly_Sugar_Data(f));
348  P.push_back(new Critical_Pair_Basic(f, StrategyFlags::SUGAR_STRATEGY));
349  }
350  // main loop
351  bool verbose = false;
352  bool very_verbose = false;
353  list<Critical_Pair_Basic *> Pnew;
354  while (not P.empty()) {
357  Critical_Pair_Basic * p = P.front();
358  report_front_pair(p, StrategyFlags::SUGAR_STRATEGY);
359  cout << "\tdegree: " << p->lcm().total_degree(0) << endl;
360  P.pop_front();
361  Pnew.push_back(p);
362  DEG_TYPE mindeg = p->lcm().total_degree();
363  for (auto pi = P.begin(); pi != P.end(); ++pi) {
364  p = *pi;
365  if (p->lcm().total_degree() < mindeg) {
366  for (auto qi = Pnew.begin(); qi != Pnew.end(); ++qi)
367  P.push_front(*qi);
368  Pnew.clear();
369  mindeg = p->lcm().total_degree();
370  report_front_pair(p, StrategyFlags::SUGAR_STRATEGY);
371  cout << "\tdegree: " << p->lcm().total_degree(0) << endl;
372  Pnew.push_back(p);
373  auto qi = pi;
374  ++qi;
375  P.erase(pi);
376  pi = qi;
377  }
378  else if (p->lcm().total_degree() == mindeg) {
379  report_front_pair(p, StrategyFlags::SUGAR_STRATEGY);
380  Pnew.push_back(p);
381  P.erase(pi);
382  }
383  }
384  // make s-poly
385  F4_Reduction_Data s(Pnew, G);
386  number_of_spolys += Pnew.size();
387  Pnew.clear();
388  if (not s.is_zero()) {
389  s.reduce_by_old();
390  s.reduce_by_new();
391  }
392  if (s.is_zero()) {
393  cout << "\tmatrix reduced to zero\n";
394  // delete s;
395  } else {
396  vector<Constant_Polynomial *> R = s.finalize();
397  for (auto r : R) {
398  cout << "\tadded " << r->leading_monomial() << endl;
399  very_verbose = false;
400  if (very_verbose) { cout << "\tadded "; r->println(); }
401  gm_update(P, G, r, StrategyFlags::SUGAR_STRATEGY);
402  }
403  }
404  }
405  cout << number_of_spolys << " s-polynomials computed and reduced\n";
406  // cleanup
407  cout << G.size() << " polynomials before interreduction\n";
408  //check_correctness(G, strategy);
409  G = reduce_basis(G);
410  cout << G.size() << " polynomials after interreduction\n";
411  list<Constant_Polynomial *> B;
412  for (Abstract_Polynomial * g : G)
413  B.push_back(new Constant_Polynomial(*g));
414  return B;
415 }
416 
417 #endif
DEG_TYPE total_degree(NVAR_TYPE m=0) const
Sum of exponents of the first m variables.
Definition: monomial.cpp:190
The general class of a polynomial.
Definition: polynomial.hpp:101
unsigned rows
number of rows in the polynomial
vector< Monomial * > M
monomials for each column
A Constant_Polynomial is a polynomial that should not change.
bool is_zero()
returns true iff all the entries are 0
vector< Abstract_Polynomial * > R
finalized list of indices of reducers for the corresponding monomials of f
void sort_pairs_by_strategy(list< T *> &P)
Applies the strategy to find the “smallest” critical pair.
virtual Polynomial_Iterator * new_iterator() const =0
An iterator that poses no risk of modifying the polynomial.
virtual void moveRight()=0
Moves right in the polynomial, to the next smaller monomial.
bool is_zero() const
Is this the additive identity?
Definition: fields.hpp:180
virtual const Prime_Field_Element & currCoeff() const =0
Reports the coefficient at the current position.
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< Monomial * > M_build
monomials while building
void print_matrix(bool show_data=false)
prints the matrix
Information necessary for a field modulo a prime.
Definition: fields.hpp:49
void set_monomial_ordering(const Monomial_Ordering *mord)
sets the Monomial_Ordering associated with this Monomial
Definition: monomial.cpp:211
virtual Monomial & leading_monomial() const =0
leading monomial – call after sort_by_order()!
vector< Poly_Sugar_Data * > strategies
strategy data for each polynomial
const Monomial & lcm() const
lcm of leading monomials of polynomials
void advance_head(unsigned i)
advances head; useful after a computation
list< Abstract_Polynomial * > R_build
indices of reducers for the corresponding elements of M
Prime_Field & ground_field()
ground field – all coefficients should be in this field
Definition: polynomial.cpp:27
void reduce_by_old()
reduces polynomials
virtual void at_generation_tasks()
hook called while first generating polynomial
Definition: strategies.hpp:88
~F4_Reduction_Data()
releases space for the matrix and deletes any strategies not already set to nullptr ...
list< Abstract_Polynomial * > reduce_basis(list< Abstract_Polynomial *>G)
Remove redundant polynomials from G.
vector< Constant_Polynomial * > finalize()
converts this to a Constant_Polynomial and returns the result
NVAR_TYPE num_vars() const
number of variables
Definition: monomial.hpp:130
virtual bool fellOff() const =0
unsigned cols
number of columns in the polynomial
Implementation of monomials.
Definition: monomial.hpp:69
Element of a field of prime characteristic.
Definition: fields.hpp:137
F4_Reduction_Data(const list< Critical_Pair_Basic *> &P, const list< Abstract_Polynomial *> &B)
encapsulation of one step of the F4 algorithm for the polynomials indicated by P and B ...
polynomial-related data for a sugar strategy
const Monomial_Ordering * mord
how the monomials are ordered
void initialize_exponents(NVAR_TYPE number_of_vars)
allocates memory for exponents This is useful when you want to allocate an array of monomials...
Definition: monomial.cpp:65
void set_strategy(Poly_Strategy_Data *psd)
sets the polynomial’s strategy to psd
Definition: polynomial.cpp:36
void add_monomials(list< Monomial *>::iterator &ti, list< Abstract_Polynomial *>::iterator &ri, const Abstract_Polynomial *g, const Monomial &u, bool new_row=false)
adds monomials of to M_build
Used to iterate through a polynomial.
Definition: polynomial.hpp:224
void initialize_many(const list< Critical_Pair_Basic *> &P)
creates the matrix
COEF_TYPE inverse() const
Returns the multiplicative inverse of this.
Definition: fields.cpp:84
virtual const Monomial & currMonomial() const =0
Reports the monomial at the current position.
void list_reducers()
lists the reducers selected for each column, in order
list< Constant_Polynomial * > f4_control(const list< Abstract_Polynomial *> &F)
equivalent to buchberger(), but for Faugère’s F4 algorithm
void reduce_by_new()
reduces polynomials
Implementation of Faugère’s F4 algorithm.
Controls the creation of s-polynomials.
Prime_Field_Element * A
coefficient data
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.
Polynomial_Ring & Rx
polynomial ring
vector< unsigned > nonzero_entries
number of nonzero entries of each row of A
vector< unsigned > heads
head of the polynomial (location of leading term)
void common_initialization(const Monomial_Ordering *ord=nullptr)
things all Monomial initializers must do
Definition: monomial.hpp:74
const list< Abstract_Polynomial * > & G
current basis of ideal
virtual Poly_Strategy_Data * strategy() const
strategy related information
Definition: polynomial.hpp:144