diff --git a/doc/source/fmpz_mod_mpoly.rst b/doc/source/fmpz_mod_mpoly.rst
index 7a40c58829..6ec3a48aab 100644
--- a/doc/source/fmpz_mod_mpoly.rst
+++ b/doc/source/fmpz_mod_mpoly.rst
@@ -628,6 +628,126 @@ Univariate Functions
Try to set *D* to the discriminant of *Ax*.
+Vectors
+--------------------------------------------------------------------------------
+
+.. type:: fmpz_mod_mpoly_vec_struct
+
+.. type:: fmpz_mod_mpoly_vec_t
+
+ A type holding a vector of :type:`fmpz_mod_mpoly_t`.
+
+.. macro:: fmpz_mod_mpoly_vec_entry(vec, i)
+
+ Macro for accessing the entry at position *i* in *vec*.
+
+.. function:: void fmpz_mod_mpoly_vec_init(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Initializes *vec* to a vector of length *len*, setting all entries to the zero polynomial.
+
+.. function:: void fmpz_mod_mpoly_vec_clear(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Clears *vec*, freeing its allocated memory.
+
+.. function:: void fmpz_mod_mpoly_vec_print(const fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Prints *vec* to standard output.
+
+.. function:: void fmpz_mod_mpoly_vec_swap(fmpz_mod_mpoly_vec_t x, fmpz_mod_mpoly_vec_t y, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Swaps *x* and *y* efficiently.
+
+.. function:: void fmpz_mod_mpoly_vec_fit_length(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Allocates room for *len* entries in *vec*.
+
+.. function:: void fmpz_mod_mpoly_vec_set(fmpz_mod_mpoly_vec_t dest, const fmpz_mod_mpoly_vec_t src, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *dest* to a copy of *src*.
+
+.. function:: void fmpz_mod_mpoly_vec_append(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Appends *f* to the end of *vec*.
+
+.. function:: slong fmpz_mod_mpoly_vec_insert_unique(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Inserts *f* without duplication into *vec* and returns its index.
+ If this polynomial already exists, *vec* is unchanged. If this
+ polynomial does not exist in *vec*, it is appended.
+
+.. function:: void fmpz_mod_mpoly_vec_set_length(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets the length of *vec* to *len*, truncating or zero-extending
+ as needed.
+
+.. function:: void fmpz_mod_mpoly_vec_randtest_not_zero(fmpz_mod_mpoly_vec_t vec, flint_rand_t state, slong len, slong poly_len, ulong exp_bound, fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *vec* to a random vector with exactly *len* entries, all nonzero,
+ with random parameters defined by *poly_len* and *exp_bound*.
+
+
+Ideals and Gröbner bases
+-------------------------------------------------------------------------------
+
+The following methods deal with ideals in `\mathbb{Z}/n\mathbb{Z}[x_1, \dots, x_m]`.
+We use monic polynomials as normalised generators.
+
+.. function:: void fmpz_mod_mpoly_spoly(fmpz_mod_mpoly_t res, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_t g, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *res* to the *S*-polynomial of *f* and *g*, scaled by making *f* and *g* monic first.
+
+.. function:: void fmpz_mod_mpoly_reduction_monic_part(fmpz_mod_mpoly_t res, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *res* to the monic remainder of multivariate
+ division with remainder with respect to the polynomials *vec*.
+
+.. function:: int fmpz_mod_mpoly_vec_is_groebner(const fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+
+ If *F* is *NULL*, checks if *G* is a Gröbner basis. If *F* is not *NULL*,
+ checks if *G* is a Gröbner basis for *F*.
+
+.. function:: int fmpz_mod_mpoly_vec_is_autoreduced(const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Checks whether the vector *F* is autoreduced (or inter-reduced).
+
+.. function:: void fmpz_mod_mpoly_vec_autoreduction(fmpz_mod_mpoly_vec_t H, const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *H* to the autoreduction (inter-reduction) of *F*.
+
+.. function:: void fmpz_mod_mpoly_vec_autoreduction_groebner(fmpz_mod_mpoly_vec_t H, const fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *H* to the autoreduction (inter-reduction) of *G*.
+ Assumes that *G* is a Gröbner basis.
+ This produces a reduced Gröbner basis, which is unique
+ (up to the sort order of the entries in the vector).
+
+.. function:: void fmpz_mod_mpoly_buchberger_naive(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+
+ Sets *G* to a Gröbner basis for *F*, computed using
+ a naive implementation of Buchberger's algorithm.
+
+.. function:: int fmpz_mod_mpoly_buchberger_naive_with_limits(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F, slong ideal_len_limit, slong poly_len_limit, const fmpz_mod_mpoly_ctx_t ctx)
+
+ As :func:`fmpz_mod_mpoly_buchberger_naive`, but halts if during the
+ execution of Buchberger's algorithm the length of the
+ ideal basis set exceeds *ideal_len_limit*, or the length of any
+ polynomial exceeds *poly_len_limit*.
+ Returns 1 for success and 0 for failure. On failure, *G* is
+ a valid basis for *F* but it might not be a Gröbner basis.
+
+
+Converting to/from other polynomial types
+-------------------------------------------------------------------------------
+
+.. function:: void fmpz_mod_mpoly_set_fmpz_mpoly(fmpz_mod_mpoly_t A, const fmpz_mpoly_t B, const fmpz_mod_mpoly_ctx_t ctxm, const fmpz_mpoly_ctx_t ctx)
+
+ Sets :type:`fmpz_mod_mpoly_t` *A* to :type:`fmpz_mpoly_t` *B* with coefficients modulo the modulus in *ctxm*.
+
+.. function:: void fmpz_mod_mpoly_get_fmpz_mpoly(fmpz_mpoly_t A, const fmpz_mod_mpoly_t B, const fmpz_mpoly_ctx_t ctx)
+
+ Sets :type:`fmpz_mpoly_t` *A* to the :type:`fmpz_mod_mpoly_t` *B*.
+
+
Internal Functions
--------------------------------------------------------------------------------
diff --git a/src/fmpz_mod_mpoly.h b/src/fmpz_mod_mpoly.h
index e103095110..f6fc6b58be 100644
--- a/src/fmpz_mod_mpoly.h
+++ b/src/fmpz_mod_mpoly.h
@@ -678,6 +678,7 @@ void fmpz_mod_mpoly_divrem_ideal_monagan_pearce(
const fmpz_mod_mpoly_t A, fmpz_mod_mpoly_struct * const * B, slong len,
const fmpz_mod_mpoly_ctx_t ctx);
+
/* Square root ***************************************************************/
int fmpz_mod_mpoly_sqrt_heap(fmpz_mod_mpoly_t Q,
@@ -711,6 +712,9 @@ int fmpz_mod_mpoly_quadratic_root(fmpz_mod_mpoly_t Q,
void fmpz_mod_mpoly_term_content(fmpz_mod_mpoly_t M,
const fmpz_mod_mpoly_t A, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_monic_part(fmpz_mod_mpoly_t res,
+ const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx);
+
int fmpz_mod_mpoly_content_vars(fmpz_mod_mpoly_t g,
const fmpz_mod_mpoly_t A, slong * vars, slong vars_length,
const fmpz_mod_mpoly_ctx_t ctx);
@@ -947,6 +951,67 @@ void fmpz_mod_mpoly_from_mpolyl_perm_inflate(fmpz_mod_mpoly_t A,
const slong * perm, const ulong * shift, const ulong * stride);
+/* Vectors of multivariate polynomials */
+
+typedef struct
+{
+ fmpz_mod_mpoly_struct * p;
+ slong alloc;
+ slong length;
+}
+fmpz_mod_mpoly_vec_struct;
+
+typedef fmpz_mod_mpoly_vec_struct fmpz_mod_mpoly_vec_t[1];
+
+#define fmpz_mod_mpoly_vec_entry(vec, i) ((vec)->p + (i))
+
+void fmpz_mod_mpoly_vec_init(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_print(const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_swap(fmpz_mod_mpoly_vec_t x, fmpz_mod_mpoly_vec_t y, const fmpz_mod_mpoly_ctx_t FLINT_UNUSED(ctx));
+void fmpz_mod_mpoly_vec_fit_length(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_clear(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_set(fmpz_mod_mpoly_vec_t dest, const fmpz_mod_mpoly_vec_t src, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_append(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx);
+slong fmpz_mod_mpoly_vec_insert_unique(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_set_length(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx);
+void fmpz_mod_mpoly_vec_randtest_not_zero(fmpz_mod_mpoly_vec_t vec, flint_rand_t state, slong len, slong poly_len, ulong exp_bound, fmpz_mod_mpoly_ctx_t ctx);
+
+
+/* Ideals and Groenber bases */
+
+void fmpz_mod_mpoly_spoly(fmpz_mod_mpoly_t res,
+ const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_t g, const fmpz_mod_mpoly_ctx_t ctx);
+
+void fmpz_mod_mpoly_reduction_monic_part(fmpz_mod_mpoly_t res,
+ const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_vec_t Iv, const fmpz_mod_mpoly_ctx_t ctx);
+
+int fmpz_mod_mpoly_vec_is_groebner(const fmpz_mod_mpoly_vec_t G,
+ const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx);
+
+void fmpz_mod_mpoly_vec_set_monic_unique(fmpz_mod_mpoly_vec_t G,
+ const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx);
+
+int fmpz_mod_mpoly_buchberger_naive_with_limits(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F,
+ slong ideal_len_limit, slong poly_len_limit, const fmpz_mod_mpoly_ctx_t ctx);
+
+void fmpz_mod_mpoly_buchberger_naive(fmpz_mod_mpoly_vec_t G,
+ const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx);
+
+int fmpz_mod_mpoly_vec_is_autoreduced(const fmpz_mod_mpoly_vec_t G,
+ const fmpz_mod_mpoly_ctx_t ctx);
+
+void fmpz_mod_mpoly_vec_autoreduction_groebner(fmpz_mod_mpoly_vec_t H,
+ const fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_ctx_t ctx);
+
+/* Convert to/from fmpz_mpoly */
+
+void fmpz_mod_mpoly_set_fmpz_mpoly(fmpz_mod_mpoly_t res,
+ const fmpz_mpoly_t f, fmpz_mod_mpoly_ctx_t ctxm, fmpz_mpoly_ctx_t ctx);
+
+void fmpz_mod_mpoly_get_fmpz_mpoly(fmpz_mpoly_t res,
+ const fmpz_mod_mpoly_t f, fmpz_mpoly_ctx_t ctx);
+
+
/******************************************************************************
Internal consistency checks
diff --git a/src/fmpz_mod_mpoly/buchberger_naive.c b/src/fmpz_mod_mpoly/buchberger_naive.c
new file mode 100644
index 0000000000..161be9efcc
--- /dev/null
+++ b/src/fmpz_mod_mpoly/buchberger_naive.c
@@ -0,0 +1,264 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz_vec.h"
+#include "fmpz_mod_mpoly.h"
+
+/* Index pairs (for Buchberger algorithm) */
+
+typedef struct
+{
+ slong a;
+ slong b;
+} pair_t;
+
+typedef struct
+{
+ pair_t * pairs;
+ slong length;
+ slong alloc;
+}
+pairs_struct;
+
+typedef pairs_struct pairs_t[1];
+
+static void
+pairs_init(pairs_t vec)
+{
+ vec->pairs = NULL;
+ vec->length = 0;
+ vec->alloc = 0;
+}
+
+static void
+pairs_fit_length(pairs_t vec, slong len)
+{
+ if (len > vec->alloc)
+ {
+ if (len < 2 * vec->alloc)
+ len = 2 * vec->alloc;
+
+ vec->pairs = flint_realloc(vec->pairs, len * sizeof(pair_t));
+ vec->alloc = len;
+ }
+}
+
+static void
+pairs_clear(pairs_t vec)
+{
+ flint_free(vec->pairs);
+}
+
+static void
+pairs_append(pairs_t vec, slong i, slong j)
+{
+ pairs_fit_length(vec, vec->length + 1);
+ vec->pairs[vec->length].a = i;
+ vec->pairs[vec->length].b = j;
+ vec->length++;
+}
+
+static pair_t
+fmpz_mod_mpoly_select_pop_pair(pairs_t pairs, const fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong len, choice, nvars;
+ pair_t result;
+
+ nvars = ctx->minfo->nvars;
+ len = pairs->length;
+ choice = 0;
+
+ if (len > 1)
+ {
+ slong i, j, a, b;
+ ulong * exp;
+ ulong * tmp;
+ ulong * lcm;
+ ulong * best_lcm;
+ ulong l, total;
+ int best;
+
+ exp = flint_malloc(sizeof(ulong) * G->length * nvars);
+ lcm = flint_malloc(sizeof(ulong) * (nvars + 1));
+ tmp = flint_malloc(sizeof(ulong) * (nvars + 1));
+ best_lcm = flint_malloc(sizeof(ulong) * (nvars + 1));
+
+ for (i = 0; i <= nvars; i++)
+ best_lcm[i] = UWORD_MAX;
+
+ for (i = 0; i < G->length; i++)
+ fmpz_mod_mpoly_get_term_exp_ui(exp + i * nvars, G->p + i, 0, ctx);
+
+ for (i = 0; i < len; i++)
+ {
+ a = pairs->pairs[i].a;
+ b = pairs->pairs[i].b;
+ total = 0;
+ best = 1;
+
+ if (ctx->minfo->ord == ORD_LEX)
+ {
+ for (j = 0; j < nvars; j++)
+ {
+ l = FLINT_MAX(exp[a * nvars + j], exp[b * nvars + j]);
+
+ if (l > best_lcm[j])
+ {
+ best = 0;
+ break;
+ }
+
+ lcm[j] = l;
+ total += l; /* total degree */
+ }
+ }
+ else /* todo: correct order */
+ {
+ for (j = 0; j < nvars; j++)
+ {
+ l = FLINT_MAX(exp[a * nvars + j], exp[b * nvars + j]);
+ total += l; /* total degree */
+
+ if (total >= best_lcm[j])
+ {
+ best = 0;
+ break;
+ }
+
+ lcm[j] = l;
+ }
+ }
+
+ if (best)
+ {
+ for (j = 0; j < nvars; j++)
+ best_lcm[j] = lcm[j];
+
+ best_lcm[nvars] = total;
+ choice = i;
+ }
+ }
+
+ flint_free(exp);
+ flint_free(tmp);
+ flint_free(lcm);
+ flint_free(best_lcm);
+ }
+
+ result = pairs->pairs[choice];
+ pairs->pairs[choice] = pairs->pairs[pairs->length - 1];
+ pairs->length--;
+
+ return result;
+}
+
+static int
+within_limits(const fmpz_mod_mpoly_t poly, slong poly_len_limit, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ if (fmpz_mod_mpoly_length(poly, ctx) > poly_len_limit)
+ return 0;
+
+ return 1;
+}
+
+static int
+fmpz_mod_mpoly_disjoint_lt(const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_t g, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ int result;
+ slong i, nvars;
+ ulong * exp1, * exp2;
+
+ nvars = ctx->minfo->nvars;
+ exp1 = flint_malloc(2 * nvars * sizeof(ulong));
+ exp2 = exp1 + nvars;
+
+ fmpz_mod_mpoly_get_term_exp_ui(exp1, f, 0, ctx);
+ fmpz_mod_mpoly_get_term_exp_ui(exp2, g, 0, ctx);
+
+ result = 1;
+ for (i = 0; i < nvars && result; i++)
+ if (exp1[i] && exp2[i])
+ result = 0;
+
+ flint_free(exp1);
+
+ return result;
+}
+
+int
+fmpz_mod_mpoly_buchberger_naive_with_limits(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F,
+ slong ideal_len_limit, slong poly_len_limit, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ pairs_t B;
+ fmpz_mod_mpoly_t h;
+ slong i, j, index_h;
+ pair_t pair;
+ int success;
+
+ fmpz_mod_mpoly_vec_set_monic_unique(G, F, ctx);
+
+ if (G->length <= 1)
+ return 1;
+
+ if (G->length >= ideal_len_limit)
+ return 0;
+
+ for (i = 0; i < G->length; i++)
+ if (!within_limits(fmpz_mod_mpoly_vec_entry(G, i), poly_len_limit, ctx))
+ return 0;
+
+ pairs_init(B);
+ fmpz_mod_mpoly_init(h, ctx);
+
+ for (i = 0; i < G->length; i++)
+ for (j = i + 1; j < G->length; j++)
+ if (!fmpz_mod_mpoly_disjoint_lt(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, j), ctx))
+ pairs_append(B, i, j);
+
+ success = 1;
+ while (B->length != 0)
+ {
+ pair = fmpz_mod_mpoly_select_pop_pair(B, G, ctx);
+
+ fmpz_mod_mpoly_spoly(h, fmpz_mod_mpoly_vec_entry(G, pair.a), fmpz_mod_mpoly_vec_entry(G, pair.b), ctx);
+ fmpz_mod_mpoly_reduction_monic_part(h, h, G, ctx);
+
+ if (!fmpz_mod_mpoly_is_zero(h, ctx))
+ {
+ /* printf("h stats %ld, %ld, %ld\n", h->length, h->bits, G->length); */
+
+ if (G->length >= ideal_len_limit || !within_limits(h, poly_len_limit, ctx))
+ {
+ success = 0;
+ break;
+ }
+
+ index_h = G->length;
+ fmpz_mod_mpoly_vec_append(G, h, ctx);
+
+ for (i = 0; i < index_h; i++)
+ if (!fmpz_mod_mpoly_disjoint_lt(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, index_h), ctx))
+ pairs_append(B, i, index_h);
+ }
+ }
+
+ fmpz_mod_mpoly_clear(h, ctx);
+ pairs_clear(B);
+
+ return success;
+}
+
+void
+fmpz_mod_mpoly_buchberger_naive(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ fmpz_mod_mpoly_buchberger_naive_with_limits(G, F, WORD_MAX, WORD_MAX, ctx);
+}
diff --git a/src/fmpz_mod_mpoly/monic_part.c b/src/fmpz_mod_mpoly/monic_part.c
new file mode 100644
index 0000000000..f91fa47536
--- /dev/null
+++ b/src/fmpz_mod_mpoly/monic_part.c
@@ -0,0 +1,42 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz.h"
+#include "fmpz_mod.h"
+#include "fmpz_vec.h"
+#include "fmpz_mod_mpoly.h"
+
+void
+fmpz_mod_mpoly_monic_part(fmpz_mod_mpoly_t res, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ if (res != f)
+ fmpz_mod_mpoly_set(res, f, ctx);
+
+ if (fmpz_mod_mpoly_is_zero(res, ctx))
+ return;
+
+ if (fmpz_sgn(res->coeffs) < 0)
+ fmpz_mod_mpoly_neg(res, res, ctx);
+
+ if (!fmpz_is_one(res->coeffs))
+ {
+ fmpz_t c;
+ fmpz_init(c);
+
+ _fmpz_vec_content(c, res->coeffs, res->length);
+ if (!fmpz_is_one(c))
+ {
+ fmpz_mod_inv(c, c, ctx->ffinfo);
+ fmpz_mod_mpoly_scalar_mul_fmpz(res, res, c, ctx);
+ }
+ fmpz_clear(c);
+ }
+}
diff --git a/src/fmpz_mod_mpoly/reduction_primitive_part.c b/src/fmpz_mod_mpoly/reduction_primitive_part.c
new file mode 100644
index 0000000000..36cf7a749a
--- /dev/null
+++ b/src/fmpz_mod_mpoly/reduction_primitive_part.c
@@ -0,0 +1,50 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz.h"
+#include "fmpz_mod_mpoly.h"
+
+void
+fmpz_mod_mpoly_reduction_monic_part(fmpz_mod_mpoly_t res, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_vec_t I, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ fmpz_t scale;
+ fmpz_mod_mpoly_struct ** Q, ** B;
+ slong i, len;
+
+ len = I->length;
+
+ fmpz_init(scale);
+ Q = flint_malloc(sizeof(fmpz_mod_mpoly_struct *) * len);
+ B = flint_malloc(sizeof(fmpz_mod_mpoly_struct *) * len);
+
+ for (i = 0; i < len; i++)
+ {
+ Q[i] = flint_malloc(sizeof(fmpz_mod_mpoly_struct));
+ fmpz_mod_mpoly_init(Q[i], ctx);
+ B[i] = fmpz_mod_mpoly_vec_entry(I, i);
+ }
+
+ fmpz_mod_mpoly_divrem_ideal(Q, res, f, B, len, ctx);
+
+ if (!fmpz_mod_mpoly_is_zero(res, ctx))
+ fmpz_mod_mpoly_make_monic(res, res, ctx);
+
+ for (i = 0; i < len; i++)
+ {
+ fmpz_mod_mpoly_clear(Q[i], ctx);
+ flint_free(Q[i]);
+ }
+
+ flint_free(Q);
+ flint_free(B);
+ fmpz_clear(scale);
+}
diff --git a/src/fmpz_mod_mpoly/spoly.c b/src/fmpz_mod_mpoly/spoly.c
new file mode 100644
index 0000000000..5541e39924
--- /dev/null
+++ b/src/fmpz_mod_mpoly/spoly.c
@@ -0,0 +1,83 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz.h"
+#include "fmpz_mod.h"
+#include "fmpz_mod_mpoly.h"
+
+void
+fmpz_mod_mpoly_spoly(fmpz_mod_mpoly_t res, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_t g, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong n, i;
+ ulong * exp, * expf, * expg;
+ fmpz_mod_mpoly_t T, U;
+
+ if (f->length == 0 || g->length == 0)
+ {
+ fmpz_mod_mpoly_zero(res, ctx);
+ return;
+ }
+
+ n = ctx->minfo->nvars;
+
+ exp = flint_malloc(sizeof(ulong) * n);
+ expf = flint_malloc(sizeof(ulong) * n);
+ expg = flint_malloc(sizeof(ulong) * n);
+ fmpz_mod_mpoly_init(T, ctx);
+ fmpz_mod_mpoly_init(U, ctx);
+
+ fmpz_mod_mpoly_get_term_exp_ui(expf, f, 0, ctx);
+ fmpz_mod_mpoly_get_term_exp_ui(expg, g, 0, ctx);
+
+ for (i = 0; i < n; i++)
+ exp[i] = FLINT_MAX(expf[i], expg[i]);
+
+ for (i = 0; i < n; i++)
+ {
+ expf[i] = exp[i] - expf[i];
+ expg[i] = exp[i] - expg[i];
+ }
+
+ if (fmpz_is_one(f->coeffs))
+ fmpz_mod_mpoly_set_coeff_ui_ui(T, 1, expf, ctx);
+ else
+ {
+ fmpz_t c;
+ fmpz_init(c);
+ fmpz_mod_inv(c, f->coeffs, ctx->ffinfo);
+ fmpz_mod_mpoly_set_coeff_fmpz_ui(T, c, expf, ctx);
+ fmpz_clear(c);
+ }
+
+ fmpz_mod_mpoly_mul(T, T, f, ctx);
+
+ if (fmpz_is_one(g->coeffs))
+ fmpz_mod_mpoly_set_coeff_ui_ui(U, 1, expg, ctx);
+ else
+ {
+ fmpz_t d;
+ fmpz_init(d);
+ fmpz_mod_inv(d, g->coeffs, ctx->ffinfo);
+ fmpz_mod_mpoly_set_coeff_fmpz_ui(U, d, expg, ctx);
+ fmpz_clear(d);
+ }
+
+ fmpz_mod_mpoly_mul(U, U, g, ctx);
+
+ fmpz_mod_mpoly_sub(res, T, U, ctx);
+
+ flint_free(exp);
+ flint_free(expf);
+ flint_free(expg);
+ fmpz_mod_mpoly_clear(T, ctx);
+ fmpz_mod_mpoly_clear(U, ctx);
+}
diff --git a/src/fmpz_mod_mpoly/test/main.c b/src/fmpz_mod_mpoly/test/main.c
index 7331db8fd3..1949194d86 100644
--- a/src/fmpz_mod_mpoly/test/main.c
+++ b/src/fmpz_mod_mpoly/test/main.c
@@ -14,6 +14,7 @@
#include "t-add_sub.c"
#include "t-add_sub_fmpz.c"
#include "t-add_sub_si.c"
+#include "t-buchberger_naive.c"
#include "t-cmp.c"
#include "t-compose_fmpz_mod_mpoly.c"
#include "t-degree.c"
@@ -55,6 +56,7 @@
#include "t-scalar_addmul_fmpz.c"
#include "t-scalar_mul_fmpz.c"
#include "t-sqrt.c"
+#include "t-to_from_fmpz_mpoly.c"
#include "t-total_degree.c"
#include "t-univar_resultant.c"
#include "t-used_vars.c"
@@ -66,6 +68,7 @@ test_struct tests[] =
TEST_FUNCTION(fmpz_mod_mpoly_add_sub),
TEST_FUNCTION(fmpz_mod_mpoly_add_sub_fmpz),
TEST_FUNCTION(fmpz_mod_mpoly_add_sub_si),
+ TEST_FUNCTION(fmpz_mod_mpoly_buchberger_naive),
TEST_FUNCTION(fmpz_mod_mpoly_compose),
TEST_FUNCTION(fmpz_mod_mpoly_cmp),
TEST_FUNCTION(fmpz_mod_mpoly_degree),
@@ -77,6 +80,7 @@ test_struct tests[] =
TEST_FUNCTION(fmpz_mod_mpoly_div_monagan_pearce),
TEST_FUNCTION(fmpz_mod_mpoly_divrem),
TEST_FUNCTION(fmpz_mod_mpoly_divrem_ideal_monagan_pearce),
+ TEST_FUNCTION(fmpz_mod_mpoly_get_fmpz_mpoly),
TEST_FUNCTION(fmpz_mod_mpoly_evaluate),
TEST_FUNCTION(fmpz_mod_mpoly_gcd_brown),
TEST_FUNCTION(fmpz_mod_mpoly_gcd_cofactors),
diff --git a/src/fmpz_mod_mpoly/test/t-buchberger_naive.c b/src/fmpz_mod_mpoly/test/t-buchberger_naive.c
new file mode 100644
index 0000000000..3114893d57
--- /dev/null
+++ b/src/fmpz_mod_mpoly/test/t-buchberger_naive.c
@@ -0,0 +1,88 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "test_helpers.h"
+#include "calcium.h"
+#include "fmpz_mod_mpoly.h"
+#include "mpoly.h"
+
+TEST_FUNCTION_START(fmpz_mod_mpoly_buchberger_naive, state)
+{
+ slong iter;
+
+ for (iter = 0; iter < 10000 * 0.1 * flint_test_multiplier(); iter++)
+ {
+ fmpz_mod_mpoly_ctx_t ctx;
+ fmpz_mod_mpoly_vec_t F, G, H;
+ slong nvars;
+ fmpz_t m;
+
+ fmpz_init(m);
+ fmpz_randtest_unsigned(m, state, n_randint(state, 2) ? 4 : n_randint(state, 100));
+ fmpz_nextprime(m, m, 0);
+ fmpz_mod_mpoly_ctx_init(ctx, 1 + n_randint(state, 4), ORD_LEX, m);
+ nvars = ctx->minfo->nvars;
+
+ fmpz_mod_mpoly_vec_init(F, 0, ctx);
+ fmpz_mod_mpoly_vec_init(G, 0, ctx);
+ fmpz_mod_mpoly_vec_init(H, 0, ctx);
+
+
+ // flint_printf("iter %wd %wd %d %{fmpz}\n\n", iter, nvars, ctx->minfo->ord, m);
+ // printf("--------------------------------------------------------------------\n");
+
+
+ if (nvars == 4)
+ fmpz_mod_mpoly_vec_randtest_not_zero(F, state, 1 + n_randint(state, 3), 1 + n_randint(state, 3), 1 + n_randint(state, 2), ctx);
+ else if (nvars == 3)
+ fmpz_mod_mpoly_vec_randtest_not_zero(F, state, 1 + n_randint(state, 4), 1 + n_randint(state, 4), 1 + n_randint(state, 2), ctx);
+ else
+ fmpz_mod_mpoly_vec_randtest_not_zero(F, state, 1 + n_randint(state, 5), 10+fmpz_bits(m), 1 + n_randint(state, 3), ctx);
+
+ // flint_printf("F = "); fmpz_mod_mpoly_vec_print(F, ctx); flint_printf("\n");
+
+ fmpz_mod_mpoly_buchberger_naive(G, F, ctx);
+
+ // flint_printf("G = "); fmpz_mod_mpoly_vec_print(G, ctx); flint_printf("\n");
+
+ if (!fmpz_mod_mpoly_vec_is_groebner(G, F, ctx))
+ {
+ flint_printf("FAIL\n\n");
+ mpoly_ordering_print(ctx->minfo->ord); printf("\n");
+ flint_printf("F = "); fmpz_mod_mpoly_vec_print(F, ctx); flint_printf("\n");
+ flint_printf("G = "); fmpz_mod_mpoly_vec_print(G, ctx); flint_printf("\n");
+ flint_abort();
+ }
+
+ fmpz_mod_mpoly_vec_autoreduction_groebner(H, G, ctx);
+
+ if (!fmpz_mod_mpoly_vec_is_groebner(H, F, ctx) || !fmpz_mod_mpoly_vec_is_autoreduced(H, ctx))
+ {
+ flint_printf("FAIL (reduced GB)\n\n");
+ mpoly_ordering_print(ctx->minfo->ord); printf("\n");
+ flint_printf("F = "); fmpz_mod_mpoly_vec_print(F, ctx); flint_printf("\n");
+ flint_printf("G = "); fmpz_mod_mpoly_vec_print(G, ctx); flint_printf("\n");
+ flint_printf("H = "); fmpz_mod_mpoly_vec_print(H, ctx); flint_printf("\n");
+ flint_abort();
+ }
+
+ fmpz_mod_mpoly_vec_clear(F, ctx);
+ fmpz_mod_mpoly_vec_clear(G, ctx);
+ fmpz_mod_mpoly_vec_clear(H, ctx);
+
+ fmpz_mod_mpoly_ctx_clear(ctx);
+
+ fmpz_clear(m);
+ }
+
+ TEST_FUNCTION_END(state);
+}
diff --git a/src/fmpz_mod_mpoly/test/t-to_from_fmpz_mpoly.c b/src/fmpz_mod_mpoly/test/t-to_from_fmpz_mpoly.c
new file mode 100644
index 0000000000..7f159f3579
--- /dev/null
+++ b/src/fmpz_mod_mpoly/test/t-to_from_fmpz_mpoly.c
@@ -0,0 +1,125 @@
+/*
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "test_helpers.h"
+#include "calcium.h"
+#include "fmpz_mod_mpoly.h"
+#include "fmpz_mpoly.h"
+#include "mpoly.h"
+
+TEST_FUNCTION_START(fmpz_mod_mpoly_get_fmpz_mpoly, state)
+{
+ slong iter;
+
+ for (iter = 0; iter < 10000 * 0.1 * flint_test_multiplier(); iter++)
+ {
+ fmpz_mpoly_ctx_t ctx;
+ fmpz_mod_mpoly_ctx_t ctxm;
+ fmpz_mod_mpoly_t A, B;
+ fmpz_mpoly_t C, D;
+ slong nvars, i;
+ fmpz_t m, c;
+
+ nvars = 1 + n_randint(state, 3);
+ fmpz_init(m);
+ fmpz_init(c);
+ fmpz_randtest_unsigned(m, state, n_randint(state, 2) ? 4 : n_randint(state, 90));
+ fmpz_nextprime(m, m, 0);
+ fmpz_mod_mpoly_ctx_init(ctxm, nvars, ORD_LEX, m);
+ fmpz_mpoly_ctx_init(ctx, nvars, ORD_LEX);
+
+ fmpz_mod_mpoly_init(A, ctxm);
+ fmpz_mod_mpoly_init(B, ctxm);
+ fmpz_mpoly_init(C, ctx);
+ fmpz_mpoly_init(D, ctx);
+
+ fmpz ** exp = (fmpz **) flint_malloc(ctx->minfo->nvars*sizeof(fmpz *));
+
+ for (i = 0; i < ctx->minfo->nvars; i++)
+ {
+ exp[i] = (fmpz *) flint_malloc(sizeof(fmpz));
+ fmpz_init(exp[i]);
+ }
+
+ // flint_printf("iter %wd %wd %d %{fmpz}\n\n", iter, nvars, ctx->minfo->ord, m);
+ // printf("--------------------------------------------------------------------\n");
+
+
+ fmpz_mpoly_randtest_bound(C, state, 1 + n_randint(state, 4), 1 + n_randint(state, 4), 1 + n_randint(state, 4), ctx);
+
+ fmpz_mod_mpoly_resize(B, C->length, ctxm);
+
+ for (i = 0; i< C->length; i++)
+ {
+ fmpz_mpoly_get_term_coeff_fmpz(c, C, i, ctx);
+ fmpz_mpoly_get_term_exp_fmpz(exp, C, i, ctx);
+
+ fmpz_mod_mpoly_set_term_exp_fmpz(B, i, exp, ctxm);
+ fmpz_mod_mpoly_set_term_coeff_fmpz(B, i, c, ctxm);
+ }
+
+ fmpz_mod_mpoly_combine_like_terms(B, ctxm);
+ fmpz_mod_mpoly_set_fmpz_mpoly(A, C, ctxm, ctx);
+
+ if (!fmpz_mod_mpoly_equal(A, B, ctxm))
+ {
+ flint_printf("FAIL\n\n");
+ mpoly_ordering_print(ctx->minfo->ord); printf("\n");
+ flint_printf("C = "); fmpz_mpoly_print_pretty(C, NULL, ctx); flint_printf("\n");
+ flint_printf("A = "); fmpz_mod_mpoly_print_pretty(A, NULL, ctxm); flint_printf("\n");
+ flint_printf("B = "); fmpz_mod_mpoly_print_pretty(B, NULL, ctxm); flint_printf("\n");
+ flint_abort();
+ }
+
+ fmpz_mpoly_resize(C, A->length, ctx);
+
+ for (i = 0; i< A->length; i++)
+ {
+ fmpz_mod_mpoly_get_term_coeff_fmpz(c, A, i, ctxm);
+ fmpz_mod_mpoly_get_term_exp_fmpz(exp, A, i, ctxm);
+
+ fmpz_mpoly_set_term_exp_fmpz(C, i, exp, ctx);
+ fmpz_mpoly_set_term_coeff_fmpz(C, i, c, ctx);
+ }
+
+ fmpz_mod_mpoly_get_fmpz_mpoly(D, A, ctx);
+
+ if (!fmpz_mpoly_equal(C, D, ctx))
+ {
+ flint_printf("FAIL\n\n");
+ mpoly_ordering_print(ctx->minfo->ord); printf("\n");
+ flint_printf("C = "); fmpz_mpoly_print_pretty(C, NULL, ctx); flint_printf("\n");
+ flint_printf("D = "); fmpz_mpoly_print_pretty(D, NULL, ctx); flint_printf("\n");
+ flint_printf("A = "); fmpz_mod_mpoly_print_pretty(A, NULL, ctxm); flint_printf("\n");
+ flint_abort();
+ }
+
+ for (i = 0; i < ctx->minfo->nvars; i++)
+ {
+ fmpz_clear(exp[i]);
+ flint_free(exp[i]);
+ }
+ flint_free(exp);
+
+ fmpz_mod_mpoly_clear(A, ctxm);
+ fmpz_mod_mpoly_clear(B, ctxm);
+ fmpz_mpoly_clear(C, ctx);
+ fmpz_mpoly_clear(D, ctx);
+
+ fmpz_mod_mpoly_ctx_clear(ctxm);
+ fmpz_mpoly_ctx_clear(ctx);
+
+ fmpz_clear(m);
+ fmpz_clear(c);
+ }
+
+ TEST_FUNCTION_END(state);
+}
diff --git a/src/fmpz_mod_mpoly/to_from_fmpz_mpoly.c b/src/fmpz_mod_mpoly/to_from_fmpz_mpoly.c
new file mode 100644
index 0000000000..2e1aa3fa5b
--- /dev/null
+++ b/src/fmpz_mod_mpoly/to_from_fmpz_mpoly.c
@@ -0,0 +1,60 @@
+/*
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+#include "fmpz.h"
+#include "fmpz_mod.h"
+#include "fmpz_mod_mpoly.h"
+#include "fmpz_mpoly.h"
+#include "mpoly.h"
+
+void
+fmpz_mod_mpoly_set_fmpz_mpoly(fmpz_mod_mpoly_t res, const fmpz_mpoly_t f, fmpz_mod_mpoly_ctx_t ctxm, fmpz_mpoly_ctx_t ctx)
+{
+ slong N = mpoly_words_per_exp(f->bits, ctx->minfo);
+ slong res_len, i;
+ FLINT_ASSERT(ctx->minfo->nvars == ctx->minfo->nvars);
+ FLINT_ASSERT(ctx->minfo->ord == ctx->minfo->ord);
+ fmpz_mod_mpoly_fit_length_reset_bits(res, f->length, f->bits, ctxm);
+ res_len = 0;
+ for (i = 0; i < f->length; i++)
+ {
+
+ fmpz_mod_set_fmpz(res->coeffs+res_len,f->coeffs + i, ctxm->ffinfo);
+
+ if (res->coeffs[res_len] == 0)
+ continue;
+
+ mpoly_monomial_set(res->exps + N*res_len, f->exps + N*i, N);
+ res_len++;
+ }
+ res->length = res_len;
+}
+
+void
+fmpz_mod_mpoly_get_fmpz_mpoly(fmpz_mpoly_t res, const fmpz_mod_mpoly_t f, fmpz_mpoly_ctx_t ctx)
+{
+ slong N = mpoly_words_per_exp(f->bits, ctx->minfo);
+ slong res_len, i;
+ FLINT_ASSERT(ctx->minfo->nvars == ctx->minfo->nvars);
+ FLINT_ASSERT(ctx->minfo->ord == ctx->minfo->ord);
+ fmpz_mpoly_fit_length_reset_bits(res, f->length, f->bits, ctx);
+ res_len = 0;
+ for (i = 0; i < f->length; i++)
+ {
+ fmpz_set(res->coeffs+res_len, f->coeffs + i);
+
+ if (res->coeffs[res_len] == 0)
+ continue;
+
+ mpoly_monomial_set(res->exps + N*res_len, f->exps + N*i, N);
+ res_len++;
+ }
+ res->length = res_len;
+}
diff --git a/src/fmpz_mod_mpoly/vec.c b/src/fmpz_mod_mpoly/vec.c
new file mode 100644
index 0000000000..62742d8a76
--- /dev/null
+++ b/src/fmpz_mod_mpoly/vec.c
@@ -0,0 +1,162 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz_mod_mpoly.h"
+
+void
+fmpz_mod_mpoly_vec_init(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ if (len == 0)
+ {
+ vec->p = NULL;
+ vec->length = 0;
+ vec->alloc = 0;
+ }
+ else
+ {
+ slong i;
+ vec->p = flint_malloc(sizeof(fmpz_mod_mpoly_struct) * len);
+ for (i = 0; i < len; i++)
+ fmpz_mod_mpoly_init(vec->p + i, ctx);
+ vec->length = vec->alloc = len;
+ }
+}
+
+void
+fmpz_mod_mpoly_vec_print(const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i;
+
+ flint_printf("[");
+ for (i = 0; i < F->length; i++)
+ {
+ fmpz_mod_mpoly_print_pretty(F->p + i, NULL, ctx);
+ if (i < F->length - 1)
+ flint_printf(", ");
+ }
+ flint_printf("]");
+}
+
+void
+fmpz_mod_mpoly_vec_swap(fmpz_mod_mpoly_vec_t x, fmpz_mod_mpoly_vec_t y, const fmpz_mod_mpoly_ctx_t FLINT_UNUSED(ctx))
+{
+ fmpz_mod_mpoly_vec_t tmp;
+ *tmp = *x;
+ *x = *y;
+ *y = *tmp;
+}
+
+void
+fmpz_mod_mpoly_vec_fit_length(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ if (len > vec->alloc)
+ {
+ slong i;
+
+ if (len < 2 * vec->alloc)
+ len = 2 * vec->alloc;
+
+ vec->p = flint_realloc(vec->p, len * sizeof(fmpz_mod_mpoly_struct));
+
+ for (i = vec->alloc; i < len; i++)
+ fmpz_mod_mpoly_init(vec->p + i, ctx);
+
+ vec->alloc = len;
+ }
+}
+
+void
+fmpz_mod_mpoly_vec_clear(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i;
+
+ for (i = 0; i < vec->alloc; i++)
+ fmpz_mod_mpoly_clear(vec->p + i, ctx);
+
+ flint_free(vec->p);
+}
+
+void
+fmpz_mod_mpoly_vec_set(fmpz_mod_mpoly_vec_t dest, const fmpz_mod_mpoly_vec_t src, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ if (dest != src)
+ {
+ slong i;
+
+ fmpz_mod_mpoly_vec_fit_length(dest, src->length, ctx);
+
+ for (i = 0; i < src->length; i++)
+ fmpz_mod_mpoly_set(dest->p + i, src->p + i, ctx);
+
+ dest->length = src->length;
+ }
+}
+
+void
+fmpz_mod_mpoly_vec_append(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ fmpz_mod_mpoly_vec_fit_length(vec, vec->length + 1, ctx);
+ fmpz_mod_mpoly_set(vec->p + vec->length, f, ctx);
+ vec->length++;
+}
+
+slong
+fmpz_mod_mpoly_vec_insert_unique(fmpz_mod_mpoly_vec_t vec, const fmpz_mod_mpoly_t f, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i;
+
+ for (i = 0; i < vec->length; i++)
+ {
+ if (fmpz_mod_mpoly_equal(vec->p + i, f, ctx))
+ return i;
+ }
+
+ fmpz_mod_mpoly_vec_append(vec, f, ctx);
+ return vec->length - 1;
+}
+
+void
+fmpz_mod_mpoly_vec_randtest_not_zero(fmpz_mod_mpoly_vec_t vec, flint_rand_t state, slong len, slong poly_len, ulong exp_bound, fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i;
+
+ fmpz_mod_mpoly_vec_set_length(vec, len, ctx);
+
+ for (i = 0; i < len; i++)
+ {
+ do {
+ fmpz_mod_mpoly_randtest_bound(vec->p + i, state, poly_len, exp_bound, ctx);
+ } while (fmpz_mod_mpoly_is_zero(vec->p + i, ctx));
+ }
+
+ vec->length = len;
+}
+
+void
+fmpz_mod_mpoly_vec_set_length(fmpz_mod_mpoly_vec_t vec, slong len, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i;
+
+ if (len > vec->length)
+ {
+ fmpz_mod_mpoly_vec_fit_length(vec, len, ctx);
+ for (i = vec->length; i < len; i++)
+ fmpz_mod_mpoly_zero(vec->p + i, ctx);
+ }
+ else if (len < vec->length)
+ {
+ for (i = len; i < vec->length; i++)
+ fmpz_mod_mpoly_zero(vec->p + i, ctx);
+ }
+
+ vec->length = len;
+}
diff --git a/src/fmpz_mod_mpoly/vec_autoreduction_groebner.c b/src/fmpz_mod_mpoly/vec_autoreduction_groebner.c
new file mode 100644
index 0000000000..357d4cbbca
--- /dev/null
+++ b/src/fmpz_mod_mpoly/vec_autoreduction_groebner.c
@@ -0,0 +1,140 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz.h"
+#include "fmpz_mod_mpoly.h"
+
+static int
+monomial_divides(const ulong * exp1, const ulong * exp2, slong nvars)
+{
+ slong i;
+ for (i = 0; i < nvars; i++)
+ if (exp1[i] > exp2[i])
+ return 0;
+ return 1;
+}
+
+void
+fmpz_mod_mpoly_vec_autoreduction_groebner(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t Gin, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i, j, nvars;
+ ulong * exp1, * exp2;
+
+ if (G != Gin)
+ fmpz_mod_mpoly_vec_set(G, Gin, ctx);
+
+ /* Content removal */
+ for (i = 0; i < G->length; i++)
+ fmpz_mod_mpoly_make_monic(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, i), ctx);
+
+ /* Make sure there are no duplicate or zero entries */
+ for (i = 0; i < G->length; i++)
+ {
+ if (fmpz_mod_mpoly_is_zero(fmpz_mod_mpoly_vec_entry(G, i), ctx))
+ {
+ fmpz_mod_mpoly_swap(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, G->length - 1), ctx);
+ fmpz_mod_mpoly_vec_set_length(G, G->length - 1, ctx);
+ }
+ else
+ {
+ for (j = i + 1; j < G->length; j++)
+ {
+ if (fmpz_mod_mpoly_equal(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, j), ctx))
+ {
+ fmpz_mod_mpoly_swap(fmpz_mod_mpoly_vec_entry(G, j), fmpz_mod_mpoly_vec_entry(G, G->length - 1), ctx);
+ fmpz_mod_mpoly_vec_set_length(G, G->length - 1, ctx);
+ }
+ }
+ }
+ }
+
+ if (G->length <= 1)
+ return;
+
+ /* First filter based on divisibility of leading terms */
+
+ nvars = ctx->minfo->nvars;
+ exp1 = flint_malloc(nvars * sizeof(ulong));
+ exp2 = flint_malloc(nvars * sizeof(ulong));
+
+ for (i = 0; i < G->length; i++)
+ {
+ fmpz_mod_mpoly_get_term_exp_ui(exp1, fmpz_mod_mpoly_vec_entry(G, i), 0, ctx);
+
+ for (j = 0; j < G->length; j++)
+ {
+ if (i != j)
+ {
+ fmpz_mod_mpoly_get_term_exp_ui(exp2, fmpz_mod_mpoly_vec_entry(G, j), 0, ctx);
+
+ if (monomial_divides(exp2, exp1, nvars))
+ {
+ fmpz_mod_mpoly_swap(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, G->length - 1), ctx);
+ fmpz_mod_mpoly_vec_set_length(G, G->length - 1, ctx);
+ break;
+ }
+ }
+ }
+ }
+
+ flint_free(exp1);
+ flint_free(exp2);
+
+ /* Now do inter-reduction */
+ if (G->length >= 2)
+ {
+ fmpz_t scale;
+ fmpz_mod_mpoly_struct ** Q, ** B;
+ slong i, j, alloc;
+
+ alloc = G->length - 1;
+
+ fmpz_init(scale);
+ Q = flint_malloc(sizeof(fmpz_mod_mpoly_struct *) * alloc);
+ B = flint_malloc(sizeof(fmpz_mod_mpoly_struct *) * alloc);
+
+ for (i = 0; i < alloc; i++)
+ {
+ Q[i] = flint_malloc(sizeof(fmpz_mod_mpoly_struct));
+ fmpz_mod_mpoly_init(Q[i], ctx);
+ }
+
+ for (i = 0; i < G->length; i++)
+ {
+ for (j = 0; j < i; j++)
+ B[j] = fmpz_mod_mpoly_vec_entry(G, j);
+ for (j = i + 1; j < G->length; j++)
+ B[j - 1] = fmpz_mod_mpoly_vec_entry(G, j);
+
+ fmpz_mod_mpoly_divrem_ideal(Q, fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, i), B, G->length - 1, ctx);
+ if (!fmpz_mod_mpoly_is_zero(fmpz_mod_mpoly_vec_entry(G, i), ctx))
+ fmpz_mod_mpoly_make_monic(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, i), ctx);
+
+ if (fmpz_mod_mpoly_is_zero(fmpz_mod_mpoly_vec_entry(G, i), ctx))
+ {
+ fmpz_mod_mpoly_swap(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, G->length - 1), ctx);
+ fmpz_mod_mpoly_vec_set_length(G, G->length - 1, ctx);
+ i--;
+ }
+ }
+
+ for (i = 0; i < alloc; i++)
+ {
+ fmpz_mod_mpoly_clear(Q[i], ctx);
+ flint_free(Q[i]);
+ }
+
+ flint_free(Q);
+ flint_free(B);
+ fmpz_clear(scale);
+ }
+}
diff --git a/src/fmpz_mod_mpoly/vec_is_autoreduced.c b/src/fmpz_mod_mpoly/vec_is_autoreduced.c
new file mode 100644
index 0000000000..369b15dda6
--- /dev/null
+++ b/src/fmpz_mod_mpoly/vec_is_autoreduced.c
@@ -0,0 +1,71 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz.h"
+#include "fmpz_mod_mpoly.h"
+
+int
+fmpz_mod_mpoly_vec_is_autoreduced(const fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i, j, len, alloc;
+ fmpz_mod_mpoly_t h;
+ fmpz_t scale;
+ fmpz_mod_mpoly_struct ** Q, ** B;
+ int result;
+
+ len = G->length;
+ alloc = len - 1;
+
+ if (len == 0)
+ return 1;
+
+ fmpz_init(scale);
+ Q = flint_malloc(sizeof(fmpz_mod_mpoly_struct *) * alloc);
+ B = flint_malloc(sizeof(fmpz_mod_mpoly_struct *) * alloc);
+
+ for (i = 0; i < alloc; i++)
+ {
+ Q[i] = flint_malloc(sizeof(fmpz_mod_mpoly_struct));
+ fmpz_mod_mpoly_init(Q[i], ctx);
+ }
+
+ fmpz_mod_mpoly_init(h, ctx);
+
+ result = 1;
+ for (i = 0; i < len && result; i++)
+ {
+ for (j = 0; j < i; j++)
+ B[j] = fmpz_mod_mpoly_vec_entry(G, j);
+ for (j = i + 1; j < G->length; j++)
+ B[j - 1] = fmpz_mod_mpoly_vec_entry(G, j);
+
+ fmpz_mod_mpoly_divrem_ideal(Q, h, fmpz_mod_mpoly_vec_entry(G, i), B, G->length - 1, ctx);
+ fmpz_mod_mpoly_make_monic(h, h, ctx);
+
+ if (fmpz_mod_mpoly_is_zero(h, ctx) || !fmpz_mod_mpoly_equal(h, fmpz_mod_mpoly_vec_entry(G, i), ctx))
+ result = 0;
+ }
+
+ for (i = 0; i < alloc; i++)
+ {
+ fmpz_mod_mpoly_clear(Q[i], ctx);
+ flint_free(Q[i]);
+ }
+
+ fmpz_mod_mpoly_clear(h, ctx);
+
+ flint_free(Q);
+ flint_free(B);
+ fmpz_clear(scale);
+
+ return result;
+}
diff --git a/src/fmpz_mod_mpoly/vec_is_groebner.c b/src/fmpz_mod_mpoly/vec_is_groebner.c
new file mode 100644
index 0000000000..54c4419b26
--- /dev/null
+++ b/src/fmpz_mod_mpoly/vec_is_groebner.c
@@ -0,0 +1,54 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz_mod_mpoly.h"
+
+int
+fmpz_mod_mpoly_vec_is_groebner(const fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i, j, len;
+ fmpz_mod_mpoly_t h;
+ int result;
+
+ len = G->length;
+
+ if (len == 0)
+ return 1;
+
+ fmpz_mod_mpoly_init(h, ctx);
+ result = 1;
+
+ for (i = 0; i < len && result; i++)
+ {
+ for (j = i + 1; j < len && result; j++)
+ {
+ fmpz_mod_mpoly_spoly(h, fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, j), ctx);
+
+ fmpz_mod_mpoly_reduction_monic_part(h, h, G, ctx);
+ if (!fmpz_mod_mpoly_is_zero(h, ctx))
+ result = 0;
+ }
+ }
+
+ if (F != NULL)
+ {
+ for (i = 0; i < F->length && result; i++)
+ {
+ fmpz_mod_mpoly_reduction_monic_part(h, fmpz_mod_mpoly_vec_entry(F, i), G, ctx);
+ if (!fmpz_mod_mpoly_is_zero(h, ctx))
+ result = 0;
+ }
+ }
+
+ fmpz_mod_mpoly_clear(h, ctx);
+ return result;
+}
diff --git a/src/fmpz_mod_mpoly/vec_set_monic_unique.c b/src/fmpz_mod_mpoly/vec_set_monic_unique.c
new file mode 100644
index 0000000000..b2b6df6494
--- /dev/null
+++ b/src/fmpz_mod_mpoly/vec_set_monic_unique.c
@@ -0,0 +1,53 @@
+/*
+ Copyright (C) 2020 Fredrik Johansson
+ Copyright (C) 2025 Andrii Yanovets
+
+ This file is part of FLINT.
+
+ FLINT is free software: you can redistribute it and/or modify it under
+ the terms of the GNU Lesser General Public License (LGPL) as published
+ by the Free Software Foundation; either version 3 of the License, or
+ (at your option) any later version. See .
+*/
+
+#include "fmpz_mod_mpoly.h"
+
+void
+fmpz_mod_mpoly_vec_set_monic_unique(fmpz_mod_mpoly_vec_t G, const fmpz_mod_mpoly_vec_t F, const fmpz_mod_mpoly_ctx_t ctx)
+{
+ slong i, j, len;
+
+ fmpz_mod_mpoly_vec_set(G, F, ctx);
+
+ len = G->length;
+
+ for (i = 0; i < len; i++)
+ {
+ /* skip zero */
+ if (fmpz_mod_mpoly_is_zero(fmpz_mod_mpoly_vec_entry(G, i), ctx))
+ {
+ fmpz_mod_mpoly_swap(fmpz_mod_mpoly_vec_entry(G, i),
+ fmpz_mod_mpoly_vec_entry(G, len - 1), ctx);
+ G->length--;
+ len--;
+ i--;
+ }
+ else
+ {
+ fmpz_mod_mpoly_make_monic(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, i), ctx);
+
+ for (j = 0; j < i; j++)
+ {
+ if (fmpz_mod_mpoly_equal(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, j), ctx))
+ {
+ fmpz_mod_mpoly_zero(fmpz_mod_mpoly_vec_entry(G, i), ctx);
+ fmpz_mod_mpoly_swap(fmpz_mod_mpoly_vec_entry(G, i), fmpz_mod_mpoly_vec_entry(G, len - 1), ctx);
+ G->length--;
+ len--;
+ i--;
+ break;
+ }
+ }
+ }
+ }
+}