// -*- C++ -*- // automatically generated by autodoc // ========== HEADER FILE src/bpol/all-irredpoly.h: ========== class all_irredpoly; // Generate all irreducible binary polynomials of given degree // ========== HEADER FILE src/bpol/bitpol-arith.h: ========== // Arithmetic with polynomials over Z/2Z // The polynomial W represented by the word w is // W = pol(w) =: \sum_k{ [bit_k(w)] * x^k} inline ulong bitpol_mult(ulong a, ulong b); // Return A * B // B=2 corresponds to multiplication with 'x' // Note that the result will silently overflow // if deg(A)+deg(B) >= BITS_PER_LONG inline ulong bitpol_square(ulong a); // Return A * A // == bit_zip(a) // == bit_zip0(a) // == bitpol_mult(a, a); inline ulong bitpol_power(ulong a, ulong e); // Return A ** e // Overflow will occur even for moderate exponents inline ulong bitpol_rem(ulong a, const ulong b); // Return R = A % B = A - (A/B)*B // Must have: B!=0 inline ulong bitpol_div(ulong a, ulong b); // Return Q = A / B // Must have B!=0. inline void bitpol_divrem(ulong a, ulong b, ulong &q, ulong &r); // Set R, Q so that A == Q * B + R. // Must have B!=0. // Equivalent to: q=bitpol_div(a,b); r=bitpol_rem(a,b); inline ulong bitpol_div_xp1(ulong a); // Return power series A / (x+1) // If A is a multiple of x+1, then the returned value // is the exact division by x+1 // Function identical to inverse_rev_gray_code(a) inline ulong bitpol_div_x2p1(ulong a); // Return power series A / (x^2+1) // If A is a multiple of x^2+1, then the returned value // is the exact division by x^2+1 // ========== HEADER FILE src/bpol/bitpol-degree.h: ========== inline ulong bitpol_deg(ulong c); // Return degree of binary polynomial C. // Return zero if C is the zero polynomial. inline ulong bitpol_h(ulong c); // Return 1-bit mask needed with modular computations. inline void bitpol_hdeg(ulong c, ulong &d, ulong &h); // Compute both degree and mask. // ========== HEADER FILE src/bpol/bitpol-deriv.h: ========== inline ulong bitpol_deriv(ulong c); // Return derivative of polynomial c // ========== HEADER FILE src/bpol/bitpol-factor.h: ========== // ----- SRCFILE=src/bpol/berlekamp.cc: ----- void setup_q_matrix(ulong c, ulong d, ulong *ss); // Compute the Q-matrix for the degree-d polynomial c. // Used in Berlekamp's factorization algorithm. ulong bitpol_refine_factors(ulong *f, ulong nf, const ulong *nn, ulong r); // Given the nullspace nn[0,...,r-1] of (Q-id) // and nf factors f[0,...,nf-1] whose product equals c // (typically nf=1 and f[0]==c) // then get all r irreducible factors of c. ulong bitpol_factor_squarefree(ulong c, ulong *f); // Fill irreducible factors of square-free polynomial c into f[] // Return number of factors. // ----- SRCFILE=src/bpol/bitpol-factor.cc: ----- ulong bitpol_factor(ulong c, ulong *f, ulong *e); // Factorize the binary polynomial c: // c = \prod_{i=0}^{fct-1}{f[i]^e[i]} // The number of factors (fct) is returned. void bitpol_sort_factorization(ulong *f, ulong *e, ulong fct); // sort wrt. to irreducible factors. ulong bitpol_test_factorization(ulong c, const ulong *f, const ulong *e, ulong fct, ulong &fi); // Test whether c == \prod_{i=0, 1, ..., fct-1}{f[i]^e[i]} // with all f[i] irreducible. // Returns zero if factorization is ok, else // the first index such that an error occured // with f[i]^e[i] is written to fi // and the return value indicates the type of error: // 1: failure with trivial factorization (c should be zero or one) // 2: factor f[fi] is reducible // 3: f[fi] does not divide c // 4: f[fi]^e[fi] does not divide c // 5: product != c // ========== HEADER FILE src/bpol/bitpol-gcd.h: ========== // GCD, EGCD and inverse with polynomials over Z/2Z // The polynomial W represented by the word w is // W = pol(w) =: \sum_k{ [bit_k(w)] * x^k} inline ulong bitpol_gcd(ulong a, ulong b); // Return polynomial gcd(A, B) inline ulong bitpol_binary_gcd(ulong a, ulong b); // Return polynomial gcd(A, B) inline ulong bitpol_egcd(ulong u, ulong v, ulong &iu, ulong &iv); // Return u3 and set u1,v1 so that gcd(u,v) == u3 == u*u1 + v*u2 // // If u3==1 then u1 is the inverse of u modulo v // and u2 is the inverse of v modulo u //. // Cf. Knuth2, p.325 // ========== HEADER FILE src/bpol/bitpol-irred.h: ========== // ----- SRCFILE=src/bpol/bitpol-irred-ben-or.cc: ----- bool bitpol_irreducible_ben_or_q(ulong c, ulong h); // Return whether C is irreducible (via the Ben-Or irreducibility test_; // h needs to be a mask with one bit set: // h == highest_one(C) >> 1 == 1UL << (degree(C)-1) // Note: routine will fail for deg(c)==BITS_PER_LONG (GCD fails) // ----- SRCFILE=src/bpol/bitpol-irred-rabin.cc: ----- bool bitpol_irreducible_rabin_q(ulong c, ulong h); // Return whether C is irreducible (via Rabin's irreducibility test). // h needs to be a mask with one bit set: // h == highest_one(C) >> 1 == 1UL << (degree(C)-1) // ----- SRCFILE=src/bpol/bitpol-spi.cc: ----- bool bitpol_need_gcd(ulong h); // Return whether GCDs are needed for irreducibility test. bool bitpol_spi_q(ulong c, ulong h); // Return whether C is a strong pseudo irreducible (SPI). // A polynomial C of degree d is an SPI if // it has no linear factors, x^(2^k)!=x for 0> 1 == 1UL << (degree(C)-1) inline bool bitpol_irreducible_q(ulong c, ulong h); inline ulong bitpol_compose_xp1(ulong c); // Return C(x+1). // Self-inverse. // If C is irreducible then C(x+1) is also irreducible. // If C is primitive then x+1 is a generator modulo C(x+1). inline ulong bitpol_recip(ulong c); // Return x^deg(C) * C(1/x) (the reciprocal polynomial) // Self-inverse. // If C is irreducible/primitive then the returned // polynomial is also irreducible/primitive. //. // Note: could use revbin(c, deg(c)-1) // ========== HEADER FILE src/bpol/bitpol-order.h: ========== // ----- SRCFILE=src/bpol/bitpol-order.cc: ----- ulong bitpol_el_order(ulong c, ulong h, const factorization &mfact, ulong a); // Return order of A modulo polynomial C. // C must be irreducible. // h must be equal 1<<(deg(c)-1). // mfact must contain the factorization of 2**deg(c)-1. // The routine may loop if either: // - the polynomial C is reducible // - h is not set correctly // - mfact is not set correctly ulong bitpol_order(ulong c, ulong h, const factorization &mfact); // Return order of polynomial C. // C must be irreducible. // h must be equal 1<<(deg(c)-1). // mfact must contain the factorization of 2**deg(c)-1. // The routine may loop if either: // - the polynomial C is reducible // - h is not set correctly // - mfact is not set correctly // Same as: bitpol_el_order(c, h, mfact, 2) // ========== HEADER FILE src/bpol/bitpol-primitive.h: ========== inline ulong test_bitpol_primitive(ulong c, ulong h, const factorization &mfact); // For an irreducible polynomial c, test whether it is primitive, // that is, whether 'x' is of maximal order. // mfact must be the factorization of the maximal order (2**n-1). // Return zero of c is primitive. // ========== HEADER FILE src/bpol/bitpol-print.h: ========== // Printing binary polynomials in different formats: // // ascii format: // bitpol_print() == bitpol_print_factorization() // x^7+x^4+x^3+1 == (x+1)^5 * (x^2+x+1) // // coefficient list: // bitpol_print_coeff() == bitpol_print_coeff_factorization() // (7,4,3,0) == (1,0)^5 * (2,1,0) // // TeX format: // bitpol_print_tex() == bitpol_print_tex_factorization() // x^{7} + x^{4} + x^{3} + 1 == \left(x + 1\right)^{5} \cdot \left(x^{2} + x + 1\right) // ----- SRCFILE=src/bpol/bitpol-print.cc: ----- void bitpol_print(const char *bla, ulong c, bool sq/*=true*/); // print as: // x^7 + x^4 + x^3 + 1 if sq==true // x^7+x^4+x^3+1 if sq==false void bitpol_print_bin(const char *bla, ulong c); // print as: 1..11.1 void bitpol_print_coeff(const char *bla, ulong c); // print as: [7, 4, 3, 0] void bitpol_print_tex(const char *bla, ulong c); // print as: x^{7} + x^{4} + x^{3} + 1 void bitpol_print_factorization(const char *bla, const ulong *f, const ulong *e, ulong fct); // print as: (x+1)^5 * (x^2+x+1) void bitpol_print_bin_factorization(const char *bla, const ulong *f, const ulong *e, ulong fct); // print as: (1.)^5 (111) void bitpol_print_coeff_factorization(const char *bla, const ulong *f, const ulong *e, ulong fct); // print as: [1,0]^5 * [2,1,0] void bitpol_print_tex_factorization(const char *bla, const ulong *f, const ulong *e, ulong fct); // print as: \left(x + 1\right)^{5} \cdot \left(x^{2} + x + 1\right) void bitpol_print_short_factorization(const char *bla, const ulong *f, const ulong *e, ulong fct); // print as: [deg1 ex1] [deg2 ex2] ... [degk exk] // ========== HEADER FILE src/bpol/bitpol-squarefree.h: ========== inline ulong bitpol_test_squarefree(ulong c); // Return 0 if polynomial is square-free // else return square factor inline ulong bitpol_pure_square_q(ulong c); // Return whether polynomial is a pure square != 1 inline ulong bitpol_pure_sqrt(ulong c); // Return t = sqrt(c) // c must be a pure square: bits at odd positions must be zero. // ----- SRCFILE=src/bpol/bitpol-squarefree.cc: ----- ulong bitpol_sreduce(ulong c); // Let c == \prod_{k}{E_k} * (\prod_{j}{S_j^{s_j}})^{2^x} // where E_k, S_j are irreducible factors // and the E_k are not necessarily distinct from the S_j. // Return polynomial f where // f = \prod_{k}{E_k} * (\prod_{j}{S_j^{s_j}})^{1} // Call repeatedly until the returned polynomial // equals c to remove all even exponents. ulong bitpol_factor_squarefree(ulong c, ulong *sf, ulong *se); // Set sf[], se[] so that // c == \prod_{i}{sf[i]^se[i]} // where all sf[i] are square-free. // Return number of factors. inline ulong bitpol_squarefree_part(ulong c); // remove all factors with even exponents inline ulong bitpol_make_squarefree(ulong c); // reduce all exponents to 1 // ========== HEADER FILE src/bpol/bitpol-srp.h: ========== inline ulong bitpol_pol2srp(ulong f, ulong d); // Return the self-reciprocal polynomial S=x^d*F(x+1/x) where d=deg(f). // W = sum(j=0, d, F(j)*x^(d-j)*(1+x^2)^j ) where // F(j) is the j-th coefficient of F. // Must have: d==degree(F) inline ulong bitpol_srp2pol(ulong s, ulong hd); // Inverse of bitpol_pol2srp(). // Must have: hd = degree(s)/2 (note: _half_ of the degree). // Only the lower half coefficients are accessed, i.e. // the routine works for degree(S) <= 2*BITS_PER_LONG-2. // ========== HEADER FILE src/bpol/bitpolmod-arith.h: ========== // Modular arithmetic with polynomials over Z/2Z. // The polynomial W represented by the word w is // W = pol(w) =: \sum_k{ [bit_k(w)] * x^k} // The modulus is C and // h needs to be a mask with one bit set: // h == highest_one(c) >> 1 == 1UL << (degree(C)-1) static inline ulong bitpolmod_times_x(ulong a, ulong c, ulong h); // Return (A * x) mod C // // If C is a primitive polynomial of degree n // successive calls will cycle through all 2**n-1 // n-bit words and the sequence of bits // (any fixed position) of a constitutes // a shift register sequence (SRS). // Start with a=2 to get an SRS that starts with // n-1 consecutive zeros (use bit 0 of a) static inline ulong bitpolmod_times_x2(ulong a, ulong c, ulong h); // Return (A * x * x) mod C // Used for squaring. static inline ulong bitpolmod_div_x(ulong a, ulong c, ulong h); // Return (A / x) mod C // C must have nonzero constant term: (c&1)==1 // // If C is a primitive polynomial of degree n // successive calls will cycle through all 2**n-1 // n-bit words and the sequence of bits // (any fixed position) of a constitutes // a shift register sequence (SRS). static inline ulong bitpolmod_inv_x(ulong c, ulong h); // Return (1 / x) mod C // C must have nonzero constant term: (c&1)==1 // c>>1 == (C-1)/x = (0-1)/x == -1/x == 1/x (mod C). static inline ulong bitpolmod_mult(ulong a, ulong b, ulong c, ulong h); // Return (A * B) mod C // // Must have deg(A) < deg(C) and deg(B) < deg(C) //. // With b=2 (== 'x') the result is identical to // bitpolmod_times_x(a, c, h) static inline ulong bitpolmod_square(ulong a, ulong c, ulong h); // Return A*A mod C // == bitpolmod_mult(a, a, c, h); // compute \sum_{k=0}^{d}{a_k\,x^k} as \sum_{k=0}^{d}{a_k\,x^{2k}} inline ulong bitpolmod_power(ulong a, ulong e, ulong c, ulong h); // Return (A ** e) mod C inline ulong bitpolmod_xpower(ulong e, ulong c, ulong h); // Return (x ** e) mod C inline ulong bitpolmod_inverse_irred(ulong a, ulong c, ulong h); // Return (A ** -1) mod C // Must have: C irreducible. // // With irreducible C the inverse of A can be obtained via // i = bitpolmod_power(a, r1, c, h) // where r1 = (h<<1)-2 = max_order - 1 = 2^degree(C) - 2 // Then 1 == bitpolmod_mult(a, i, c, h) inline ulong bitpolmod_inverse(ulong a, ulong c); // Return the inverse of A modulo C if it exists, else zero. // Must have deg(A) < deg(C) inline ulong bitpolmod_divide(ulong a, ulong b, ulong c, ulong h); // Return A/B modulo C. // Must have: gcd(b,c)==1 static inline ulong bitpolmod_sqrt(ulong a, ulong c, ulong h); // Return sqrt(A) mod C // Using the identity sqrt(A) = A^(2^(n-1)) where n=deg(C) // ========== HEADER FILE src/bpol/bitpolmod-minpoly.h: ========== // ----- SRCFILE=src/bpol/bitpolmod-minpoly.cc: ----- ulong bitpolmod_minpoly(ulong a, ulong c, ulong n, ulong &bp); // Compute the minimal polynomial p(x) of A modulo C // Return the degree of p(). // The polynomial p() is written to bp // (its coefficients are in GF(2)). ulong bitpolmod_minpoly2(ulong a, ulong c, ulong n, ulong &bp); // Compute the minimal polynomial p(x) of A modulo C // Return the degree of p(). // The polynomial p() is written to bp // (its coefficients are in GF(2)). // NOTE: prefer the routine bitpolmod_minpoly()! // ========== HEADER FILE src/bpol/bitpolmod-solvequadratic.h: ========== // ----- SRCFILE=src/bpol/bitpolmod-solvequadratic.cc: ----- bool bitpolmod_solve_reduced_quadratic(ulong c, ulong& r, ulong m); // Solve z^2 + z == c modulo m // Return whether solutions exist. // If so, one solutions is written to r. // The other solution is r+1. bool bitpolmod_solve_quadratic(ulong a, ulong b, ulong c, ulong& r0, ulong& r1, ulong m); // Solve a*z^2 + b*z + c == 0 modulo m // Return whether solutions exist. // If so, the solutions are written to r0 and r1. // Must have: m irreducible. // ========== HEADER FILE src/bpol/clhca.h: ========== // Cyclic (additive) Linear Hybrid Cellular Automata (CLHCA). inline ulong clhca_next(ulong x, ulong r, ulong n); // Compute next state of a linear hybrid cellular automaton // with cyclic boundary condition (CLHCA). inline ulong clhca2poly(ulong r, ulong n); // Compute the binary polynomial corresponding // to the length-n CLHCA with rule r. inline ulong clhca2poly_too(ulong r, ulong n); // Compute the binary polynomial corresponding // to the length-n CLHCA with rule r. // This is the more general (and expensive) version. // ========== HEADER FILE src/bpol/fcsr.h: ========== class fcsr; // Feedback Carry Shift Register (FCSR) // Produces a shift register sequence (SRS) // generated by a_k mod c where // c should be a prime that has the primitive root 2 // We have a_k = a_0 * 2^k (mod c) [by default a_0=1] // The period of the FCSR is the order of two mod c // ( ==c-1 if c is prime with primitive root 2) // ========== HEADER FILE src/bpol/gf2n-trace.h: ========== inline ulong gf2n_fast_trace(ulong a, ulong tv); // Fast computation of the trace of a, // using the pre-calculated table tv from gf2n_trace_vector() // ----- SRCFILE=src/bpol/gf2n-trace.cc: ----- ulong gf2n_trace(ulong a, ulong c, ulong h); // Return the trace of A, an element of GF(2**k). // GF(2**k) represented as polynomials modulo C. // trace(A) = \sum_{j=0}^{k-1}{A^(2^j)} // h needs to be a mask with one bit set: // h == highest_one(c) >> 1 == 1UL << (degree(C)-1) // the value returned is zero or one ulong gf2n_trace_vector(ulong g, ulong c, ulong h); // Return traces of powers of G, an element of GF(2**k). // GF(2**k) represented as polynomials modulo C. // h needs to be a mask with one bit set: // h == highest_one(c) >> 1 == 1UL << (degree(C)-1) // The returned word can be used for accelerated // trace computations if G is a generator of GF(2**k) // represented as polynomials modulo (irreducible) C // (usually G='x'=2UL and C a primitive polynomial) ulong gf2n_trace_vector_x(ulong c, ulong n); // Return vector of traces of powers of x, where // x is a root of the irreducible polynomial C. // Must have: n == degree(C) ulong gf2n_half_trace(ulong a, ulong c, ulong h); // Return the half-trace of A, an element of GF(2**k). // GF(2**k) represented as polynomials modulo C. // k must be odd. // halftrace(A) = \sum_{j=0}^{(k-1)/2}{A^(4^(j))} // Let T be the trace, H the halftrace of A. // Then: H + H^2 = A + T // therefore: // if (T==0) both H and H+1 solve y^2+y=A; // else (T==1) y^2+y=A has no solutions. // further: // H(A) + H(A)^2 == H(A + A^2) // h needs to be a mask with one bit set: // h == highest_one(c) >> 1 == 1UL << (degree(C)-1) // ========== HEADER FILE src/bpol/gf2n.h: ========== // ----- SRCFILE=src/bpol/gf2n-minpoly.cc: ----- ulong gf2n_minpoly(GF2n a, ulong &bp); // Compute the minimal polynomial p(x) of a \in GF(2**n). // Return the degree of p(). // The polynomial p() is written to bp // (its coefficients are in GF(2)). // A version that does not depend on class GF2n is // given in bpol/bitpolmod-minpoly.cc ulong gf2n_minpoly2(GF2n a, ulong &bp); // Compute the minimal polynomial p(x) of a \in GF(2**n). // Return the degree of p(). // The polynomial p() is written to bp // (its coefficients are in GF(2)). // A version that does not depend on class GF2n is // given in bpol/bitpolmod-minpoly.cc // NOTE: prefer the routine gf2n_minpoly()! GF2n gf2n_eval_poly(GF2n a, ulong bp); // Evalute polynomial bp (coefficients \in GF(2)) // for argument a \in GF(2**n). // ----- SRCFILE=src/bpol/gf2n-solvequadratic.cc: ----- bool gf2n_solve_reduced_quadratic(GF2n c, GF2n& r); // Solve z^2 + z == c // Return whether solutions exist. // If so, one solutions is written to r. // The other solution is r+1. bool gf2n_solve_quadratic(GF2n a, GF2n b, GF2n c, GF2n& r0, GF2n& r1); // Solve a*z^2 + b*z + c == 0 // Return whether solutions exist. // If so, the solutions are written to r0 and r1. // ----- SRCFILE=src/bpol/gf2n-order.cc: ----- ulong gf2n_order(ulong g, ulong c, ulong h, const factorization &mfact); // Return order of g \in GF(2**n) with field polynomial c. // c must be irreducible // h must be equal 1<<(deg(c)-1) // mfact must contain the factorization of 2**deg(c)-1 // The routine may loop if either: // - the polynomial c is reducible // - deg(g) >= deg(c) // - h is not set correctly // - mfact is not set correctly class GF2n; // Implementation of binary finite fields GF(2**n) // with the arithmetic operations. // GF2n::init(n) _MUST_ be called before usage. // n must be <= BITS_PER_LONG. // ========== HEADER FILE src/bpol/lfsr.h: ========== class lfsr; // (binary) Linear Feedback Shift Register // Produces a shift register sequence (SRS) // generated by a_k=x^k (mod c) where // c is a primitive polynomial of degree n. // The period of SRS is 2^n - 1 // (non-primitive c lead to smaller periods) // ========== HEADER FILE src/bpol/lfsr64.h: ========== class lfsr64; // (binary) Linear Feedback Shift Register // Produces a shift register sequence (SRS) // generated by a_k mod c where // c defaults to the minimum weight primitive polynomial of degree 64: // c = x^64 + x^4 + x^3 + x + 1 // and // a_k = x^k (mod c) // period of lfsr is 2^64 - 1 // ========== HEADER FILE src/bpol/lhca.h: ========== inline ulong lhca_next(ulong x, ulong r, ulong m); // LHCA := (1-dim) Linear Hybrid Cellular Automaton. // Return next state (after x) of the LHCA with // rule (defined by) r: // Rule 150 is applied for cells where r is one, rule 90 else. // Rule 150 := next(x) = x + leftbit(x) + rightbit(x) // Rule 90 := next(x) = leftbit(x) + rightbit(x) // Length defined by m: // m has to be a burst of the n lowest bits (n: length of automaton) inline ulong lhca2poly(ulong r, ulong n); // Return binary polynomial p that corresponds to the length-n LHCA rule r. // The returned polynomial p has degree n. // // If r gives maximal period then p is primitive. // // Algorithm: // Let r=[r(n-1),...,r(2),r(1),r(0)] // Initialize: p2=0, p1=1 // Iterate for k=0...n-1: {p1, p2} = { (x+r)*p1+p2 ,p1} // ----- SRCFILE=src/bpol/bitpol2lhca.cc: ----- ulong poly2lhca(ulong p); // Return LHCA rule corresponding to the binary polynomial P. // Must have: P irreducible. // ========== HEADER FILE src/bpol/mersenne-coprime.h: ========== class mersenne_coprime; // For k=1 .. m, determine whether gcd(k,m) != 0 where m=2^e-1. // Must have e <= 64 (=sizeof(umod_t)). // ========== HEADER FILE src/bpol/necklace2bitpol.h: ========== class necklace2bitpol; // Generate irreducible polynomials from necklaces. // ========== HEADER FILE src/bpol/normal-solvequadratic.h: ========== inline ulong normal_solve_reduced_quadratic(ulong c); // Solve x^2 + x = c // Must have: trace(c)==0, i.e. parity(c)==0 // Return one solution x, the other solution equals 1+x, // that is, the complement of x. inline ulong normal_solve_reduced_quadratic_q(ulong c, ulong &x); // Return t, the trace of c. // If t==0 then x^2 + x = c is solvable // and a solution is written to x. // ========== HEADER FILE src/bpol/normalbasis.h: ========== // ----- SRCFILE=src/bpol/bitpol-normal.cc: ----- bool bitpol_normal2_q(ulong c, ulong n); // Return whether polynomial c (of degree n) is normal. // Must have: c irreducible. bool bitpol_normal_q(ulong c, ulong n, ulong iq/*=nullptr*/, ulong *M/*=nullptr*/); // Return whether polynomial c (of degree n) is normal. // Set iq to 1 for polynomials known to be irreducible. // If M is given then the multiplication matrix is computed. // ----- SRCFILE=src/bpol/normal-mult.cc: ----- ulong normal_mult(ulong a, ulong b, const ulong *M, ulong n); // Multiply two elements (a and b in GF(2^n)) in normal basis representation. // The multiplication matrix has to be supplied in M. // ========== HEADER FILE src/bpol/normalpoly-dual.h: ========== // ----- SRCFILE=src/bpol/normalpoly-dual.cc: ----- ulong gf2n_xx2k_trace(ulong c, ulong deg); // Return vector T of traces T[k]=trace(ek), // where ek = x*x^(2^k), k=0..deg-1, and // x is a root of the irreducible polynomial C. // Must have: deg == degree(C) // The polynomial C is normal if and only if // gcd(x^deg-1, T)==1 where T is taken as a polynomial. // The dual basis is generated by sum(k=0..deg-1, D[k]*x^(2^k)) where // D = T^-1 mod x^deg-1. The dual basis is also normal. // If the only nonzero component of T is T[0], then // the basis is self-dual. ulong gf2n_dual_normal(ulong c, ulong deg, ulong ntc/*=0*/, ulong *ntd/*=0*/); // Return the minimal polynomial CS for the dual (normal) basis // with the irreducible normal polynomial C. // Return zero if C is not normal. // Must have: deg == degree(C). // If ntc is supplied it must be equal to gf2n_xx2k_trace(c, deg). // If ntd is nonzero, ntc^-1 (mod x^deg-1) is written to it. // ========== HEADER FILE src/bpol/num-bitpol.h: ========== // num_irredpol_tab[n] == number of degree-n binary irreducible polynomials // num_primpol_tab[n] == number of degree-n binary primitive polynomials // num_normalpol_tab[n] == number of degree-n binary normal polynomials // ========== HEADER FILE src/bpol/poly-tab.h: ==========