source: git/ntl/include/NTL/ZZ_pX.h @ 92c0ac

spielwiese
Last change on this file since 92c0ac was 92c0ac, checked in by Hans Schönemann <hannes@…>, 21 years ago
NTL-5.2 git-svn-id: file:///usr/local/Singular/svn/trunk@6320 2c84dea3-7e68-4137-9b89-c4e89433aadc
  • Property mode set to 100644
File size: 34.3 KB
Line 
1
2
3#ifndef NTL_ZZ_pX__H
4#define NTL_ZZ_pX__H
5
6#include <NTL/vector.h>
7#include <NTL/ZZ_p.h>
8#include <NTL/vec_ZZ.h>
9#include <NTL/vec_ZZ_p.h>
10#include <NTL/FFT.h>
11
12NTL_OPEN_NNS
13
14
15
16
17// some cross-over points
18// macros are used so as to be consistent with zz_pX
19
20#define NTL_ZZ_pX_FFT_CROSSOVER (20) 
21#define NTL_ZZ_pX_NEWTON_CROSSOVER (45)
22#define NTL_ZZ_pX_DIV_CROSSOVER (90)
23#define NTL_ZZ_pX_HalfGCD_CROSSOVER (25)
24#define NTL_ZZ_pX_GCD_CROSSOVER (180)
25#define NTL_ZZ_pX_BERMASS_CROSSOVER (90)
26#define NTL_ZZ_pX_TRACE_CROSSOVER (90)
27
28
29
30/************************************************************
31
32                         ZZ_pX
33
34The class ZZ_pX implements polynomial arithmetic modulo p.
35Polynomials are represented as vec_ZZ_p's.
36If f is a ZZ_pX, then f.rep is a vec_ZZ_p.
37The zero polynomial is represented as a zero length vector.
38Otherwise. f.rep[0] is the constant-term, and f.rep[f.rep.length()-1]
39is the leading coefficient, which is always non-zero.
40The member f.rep is public, so the vector representation is fully
41accessible.
42Use the member function normalize() to strip leading zeros.
43
44**************************************************************/
45
46
47class ZZ_pX {
48
49public:
50
51typedef vec_ZZ_p VectorBaseType; 
52
53
54vec_ZZ_p rep;
55
56
57/***************************************************************
58
59          Constructors, Destructors, and Assignment
60
61****************************************************************/
62
63
64ZZ_pX()
65//  initial value 0
66
67   { }
68
69
70ZZ_pX(INIT_SIZE_TYPE, long n) { rep.SetMaxLength(n); }
71
72ZZ_pX(const ZZ_pX& a) : rep(a.rep) { }
73// initial value is a
74
75
76ZZ_pX& operator=(const ZZ_pX& a) 
77   { rep = a.rep; return *this; }
78
79~ZZ_pX() { }
80
81void normalize();
82// strip leading zeros
83
84void SetMaxLength(long n) 
85// pre-allocate space for n coefficients.
86// Value is unchanged
87
88   { rep.SetMaxLength(n); }
89
90
91void kill() 
92// free space held by this polynomial.  Value becomes 0.
93
94   { rep.kill(); }
95
96static const ZZ_pX& zero();
97
98
99ZZ_pX(ZZ_pX& x, INIT_TRANS_TYPE) : rep(x.rep, INIT_TRANS) { }
100
101inline ZZ_pX(long i, const ZZ_p& c);
102inline ZZ_pX(long i, long c);
103
104ZZ_pX& operator=(long a);
105ZZ_pX& operator=(const ZZ_p& a);
106
107
108};
109
110
111
112
113/********************************************************************
114
115                           input and output
116
117I/O format:
118
119   [a_0 a_1 ... a_n],
120
121represents the polynomial a_0 + a_1*X + ... + a_n*X^n.
122
123On output, all coefficients will be integers between 0 and p-1,
124amd a_n not zero (the zero polynomial is [ ]).
125On input, the coefficients are arbitrary integers which are
126then reduced modulo p, and leading zeros stripped.
127
128*********************************************************************/
129
130
131NTL_SNS istream& operator>>(NTL_SNS istream& s, ZZ_pX& x);
132NTL_SNS ostream& operator<<(NTL_SNS ostream& s, const ZZ_pX& a);
133
134
135
136
137/**********************************************************
138
139                   Some utility routines
140
141***********************************************************/
142
143
144inline long deg(const ZZ_pX& a) { return a.rep.length() - 1; }
145// degree of a polynomial.
146// note that the zero polynomial has degree -1.
147
148const ZZ_p& coeff(const ZZ_pX& a, long i);
149// zero if i not in range
150
151void GetCoeff(ZZ_p& x, const ZZ_pX& a, long i);
152// x = a[i], or zero if i not in range
153
154const ZZ_p& LeadCoeff(const ZZ_pX& a);
155// zero if a == 0
156
157const ZZ_p& ConstTerm(const ZZ_pX& a);
158// zero if a == 0
159
160void SetCoeff(ZZ_pX& x, long i, const ZZ_p& a);
161// x[i] = a, error is raised if i < 0
162
163void SetCoeff(ZZ_pX& x, long i, long a);
164
165void SetCoeff(ZZ_pX& x, long i);
166// x[i] = 1, error is raised if i < 0
167
168inline ZZ_pX::ZZ_pX(long i, const ZZ_p& a)
169   { SetCoeff(*this, i, a); } 
170
171inline ZZ_pX::ZZ_pX(long i, long a)
172   { SetCoeff(*this, i, a); } 
173
174void SetX(ZZ_pX& x);
175// x is set to the monomial X
176
177long IsX(const ZZ_pX& a);
178// test if x = X
179
180inline void clear(ZZ_pX& x) 
181// x = 0
182
183   { x.rep.SetLength(0); }
184
185inline void set(ZZ_pX& x)
186// x = 1
187
188   { x.rep.SetLength(1); set(x.rep[0]); }
189
190inline void swap(ZZ_pX& x, ZZ_pX& y)
191// swap x & y (only pointers are swapped)
192
193   { swap(x.rep, y.rep); }
194
195void random(ZZ_pX& x, long n);
196inline ZZ_pX random_ZZ_pX(long n)
197   { ZZ_pX x; random(x, n); NTL_OPT_RETURN(ZZ_pX, x); }
198// generate a random polynomial of degree < n
199
200void trunc(ZZ_pX& x, const ZZ_pX& a, long m);
201// x = a % X^m
202
203inline ZZ_pX trunc(const ZZ_pX& a, long m)
204   { ZZ_pX x; trunc(x, a, m); NTL_OPT_RETURN(ZZ_pX, x); }
205
206void RightShift(ZZ_pX& x, const ZZ_pX& a, long n);
207// x = a/X^n
208
209inline ZZ_pX RightShift(const ZZ_pX& a, long n)
210   { ZZ_pX x; RightShift(x, a, n); NTL_OPT_RETURN(ZZ_pX, x); }
211
212void LeftShift(ZZ_pX& x, const ZZ_pX& a, long n);
213// x = a*X^n
214
215inline ZZ_pX LeftShift(const ZZ_pX& a, long n)
216   { ZZ_pX x; LeftShift(x, a, n); NTL_OPT_RETURN(ZZ_pX, x); }
217
218#ifndef NTL_TRANSITION
219
220inline ZZ_pX operator>>(const ZZ_pX& a, long n)
221   { ZZ_pX x; RightShift(x, a, n); NTL_OPT_RETURN(ZZ_pX, x); }
222
223inline ZZ_pX operator<<(const ZZ_pX& a, long n)
224   { ZZ_pX x; LeftShift(x, a, n); NTL_OPT_RETURN(ZZ_pX, x); }
225
226inline ZZ_pX& operator<<=(ZZ_pX& x, long n)
227   { LeftShift(x, x, n); return x; }
228
229inline ZZ_pX& operator>>=(ZZ_pX& x, long n)
230   { RightShift(x, x, n); return x; }
231
232#endif
233
234
235
236void diff(ZZ_pX& x, const ZZ_pX& a);
237// x = derivative of a
238
239inline ZZ_pX diff(const ZZ_pX& a)
240   { ZZ_pX x; diff(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
241
242
243void MakeMonic(ZZ_pX& x);
244
245void reverse(ZZ_pX& c, const ZZ_pX& a, long hi);
246
247inline ZZ_pX reverse(const ZZ_pX& a, long hi)
248   { ZZ_pX x; reverse(x, a, hi); NTL_OPT_RETURN(ZZ_pX, x); }
249
250inline void reverse(ZZ_pX& c, const ZZ_pX& a)
251{  reverse(c, a, deg(a)); }
252
253inline ZZ_pX reverse(const ZZ_pX& a)
254   { ZZ_pX x; reverse(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
255
256inline void VectorCopy(vec_ZZ_p& x, const ZZ_pX& a, long n)
257   { VectorCopy(x, a.rep, n); }
258
259inline vec_ZZ_p VectorCopy(const ZZ_pX& a, long n)
260   { return VectorCopy(a.rep, n); }
261
262
263
264
265/*******************************************************************
266
267                        conversion routines
268
269********************************************************************/
270
271
272
273void conv(ZZ_pX& x, long a);
274void conv(ZZ_pX& x, const ZZ& a);
275void conv(ZZ_pX& x, const ZZ_p& a);
276void conv(ZZ_pX& x, const vec_ZZ_p& a);
277
278inline ZZ_pX to_ZZ_pX(long a)
279   { ZZ_pX x; conv(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
280
281inline ZZ_pX to_ZZ_pX(const ZZ& a)
282   { ZZ_pX x; conv(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
283
284inline ZZ_pX to_ZZ_pX(const ZZ_p& a)
285   { ZZ_pX x; conv(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
286
287inline ZZ_pX to_ZZ_pX(const vec_ZZ_p& a)
288   { ZZ_pX x; conv(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
289
290
291
292/*************************************************************
293
294                        Comparison
295
296**************************************************************/
297
298long IsZero(const ZZ_pX& a); 
299
300long IsOne(const ZZ_pX& a);
301
302inline long operator==(const ZZ_pX& a, const ZZ_pX& b)
303{
304   return a.rep == b.rep;
305}
306
307inline long operator!=(const ZZ_pX& a, const ZZ_pX& b)
308{
309   return !(a == b);
310}
311
312long operator==(const ZZ_pX& a, long b);
313long operator==(const ZZ_pX& a, const ZZ_p& b);
314
315inline long operator==(long a, const ZZ_pX& b) { return b == a; }
316inline long operator==(const ZZ_p& a, const ZZ_pX& b) { return b == a; }
317
318inline long operator!=(const ZZ_pX& a, long b) { return !(a == b); }
319inline long operator!=(const ZZ_pX& a, const ZZ_p& b) { return !(a == b); }
320
321inline long operator!=(long a, const ZZ_pX& b) { return !(a == b); }
322inline long operator!=(const ZZ_p& a, const ZZ_pX& b) { return !(a == b); }
323
324
325/***************************************************************
326
327                         Addition
328
329****************************************************************/
330
331void add(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
332// x = a + b
333
334void sub(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
335// x = a - b
336
337void negate(ZZ_pX& x, const ZZ_pX& a);
338// x = -a
339
340// scalar versions
341
342void add(ZZ_pX& x, const ZZ_pX& a, const ZZ_p& b); // x = a + b
343void add(ZZ_pX& x, const ZZ_pX& a, long b);
344
345inline void add(ZZ_pX& x, const ZZ_p& a, const ZZ_pX& b) { add(x, b, a); }
346inline void add(ZZ_pX& x, long a, const ZZ_pX& b) { add(x, b, a); }
347
348
349void sub(ZZ_pX & x, const ZZ_pX& a, const ZZ_p& b); // x = a - b
350
351void sub(ZZ_pX& x, const ZZ_pX& a, long b);
352void sub(ZZ_pX& x, const ZZ_pX& a, const ZZ_p& b);
353
354void sub(ZZ_pX& x, long a, const ZZ_pX& b);
355void sub(ZZ_pX& x, const ZZ_p& a, const ZZ_pX& b);
356
357inline ZZ_pX operator+(const ZZ_pX& a, const ZZ_pX& b)
358   { ZZ_pX x; add(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
359
360inline ZZ_pX operator+(const ZZ_pX& a, const ZZ_p& b)
361   { ZZ_pX x; add(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
362
363inline ZZ_pX operator+(const ZZ_pX& a, long b)
364   { ZZ_pX x; add(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
365
366inline ZZ_pX operator+(const ZZ_p& a, const ZZ_pX& b)
367   { ZZ_pX x; add(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
368
369inline ZZ_pX operator+(long a, const ZZ_pX& b)
370   { ZZ_pX x; add(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
371
372
373inline ZZ_pX operator-(const ZZ_pX& a, const ZZ_pX& b)
374   { ZZ_pX x; sub(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
375
376inline ZZ_pX operator-(const ZZ_pX& a, const ZZ_p& b)
377   { ZZ_pX x; sub(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
378
379inline ZZ_pX operator-(const ZZ_pX& a, long b)
380   { ZZ_pX x; sub(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
381
382inline ZZ_pX operator-(const ZZ_p& a, const ZZ_pX& b)
383   { ZZ_pX x; sub(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
384
385inline ZZ_pX operator-(long a, const ZZ_pX& b)
386   { ZZ_pX x; sub(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
387
388
389inline ZZ_pX& operator+=(ZZ_pX& x, const ZZ_pX& b)
390   { add(x, x, b); return x; }
391
392inline ZZ_pX& operator+=(ZZ_pX& x, const ZZ_p& b)
393   { add(x, x, b); return x; }
394
395inline ZZ_pX& operator+=(ZZ_pX& x, long b)
396   { add(x, x, b); return x; }
397
398inline ZZ_pX& operator-=(ZZ_pX& x, const ZZ_pX& b)
399   { sub(x, x, b); return x; }
400
401inline ZZ_pX& operator-=(ZZ_pX& x, const ZZ_p& b)
402   { sub(x, x, b); return x; }
403
404inline ZZ_pX& operator-=(ZZ_pX& x, long b)
405   { sub(x, x, b); return x; }
406
407
408inline ZZ_pX operator-(const ZZ_pX& a) 
409   { ZZ_pX x; negate(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
410
411inline ZZ_pX& operator++(ZZ_pX& x) { add(x, x, 1); return x; }
412inline void operator++(ZZ_pX& x, int) { add(x, x, 1); }
413inline ZZ_pX& operator--(ZZ_pX& x) { sub(x, x, 1); return x; }
414inline void operator--(ZZ_pX& x, int) { sub(x, x, 1); }
415
416/*****************************************************************
417
418                        Multiplication
419
420******************************************************************/
421
422
423void mul(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
424// x = a * b
425
426void sqr(ZZ_pX& x, const ZZ_pX& a);
427inline ZZ_pX sqr(const ZZ_pX& a)
428   { ZZ_pX x; sqr(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
429// x = a^2
430
431void mul(ZZ_pX & x, const ZZ_pX& a, const ZZ_p& b); 
432void mul(ZZ_pX& x, const ZZ_pX& a, long b);
433
434
435inline void mul(ZZ_pX& x, const ZZ_p& a, const ZZ_pX& b)
436   { mul(x, b, a); }
437
438inline void mul(ZZ_pX& x, long a, const ZZ_pX& b)
439   { mul(x, b, a); }
440
441inline ZZ_pX operator*(const ZZ_pX& a, const ZZ_pX& b)
442   { ZZ_pX x; mul(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
443
444inline ZZ_pX operator*(const ZZ_pX& a, const ZZ_p& b)
445   { ZZ_pX x; mul(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
446
447inline ZZ_pX operator*(const ZZ_pX& a, long b)
448   { ZZ_pX x; mul(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
449
450inline ZZ_pX operator*(const ZZ_p& a, const ZZ_pX& b)
451   { ZZ_pX x; mul(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
452
453inline ZZ_pX operator*(long a, const ZZ_pX& b)
454   { ZZ_pX x; mul(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
455
456inline ZZ_pX& operator*=(ZZ_pX& x, const ZZ_pX& b)
457   { mul(x, x, b); return x; }
458
459inline ZZ_pX& operator*=(ZZ_pX& x, const ZZ_p& b)
460   { mul(x, x, b); return x; }
461
462inline ZZ_pX& operator*=(ZZ_pX& x, long b)
463   { mul(x, x, b); return x; }
464
465
466void PlainMul(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
467// always uses the "classical" algorithm
468
469void PlainSqr(ZZ_pX& x, const ZZ_pX& a);
470// always uses the "classical" algorithm
471
472
473void FFTMul(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
474// always uses the FFT
475
476void FFTSqr(ZZ_pX& x, const ZZ_pX& a);
477// always uses the FFT
478
479void MulTrunc(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b, long n);
480// x = a * b % X^n
481
482inline ZZ_pX MulTrunc(const ZZ_pX& a, const ZZ_pX& b, long n)
483   { ZZ_pX x; MulTrunc(x, a, b, n); NTL_OPT_RETURN(ZZ_pX, x); }
484
485void PlainMulTrunc(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b, long n);
486void FFTMulTrunc(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b, long n);
487
488void SqrTrunc(ZZ_pX& x, const ZZ_pX& a, long n);
489// x = a^2 % X^n
490
491inline ZZ_pX SqrTrunc(const ZZ_pX& a, long n)
492   { ZZ_pX x; SqrTrunc(x, a, n); NTL_OPT_RETURN(ZZ_pX, x); }
493
494void PlainSqrTrunc(ZZ_pX& x, const ZZ_pX& a, long n);
495void FFTSqrTrunc(ZZ_pX& x, const ZZ_pX& a, long n);
496
497
498void power(ZZ_pX& x, const ZZ_pX& a, long e);
499inline ZZ_pX power(const ZZ_pX& a, long e)
500   { ZZ_pX x; power(x, a, e); NTL_OPT_RETURN(ZZ_pX, x); }
501
502
503// The following data structures and routines allow one
504// to hand-craft various algorithms, using the FFT convolution
505// algorithms directly.
506// Look in the file ZZ_pX.c for examples.
507
508
509
510
511// FFT representation of polynomials
512
513class FFTRep {
514public:
515   long k;                // a 2^k point representation
516   long MaxK;             // maximum space allocated
517   long **tbl;
518   long NumPrimes; 
519
520   void SetSize(long NewK);
521
522   FFTRep(const FFTRep& R);
523   FFTRep& operator=(const FFTRep& R);
524
525   FFTRep() { k = MaxK = -1; tbl = 0; NumPrimes = 0; }
526   FFTRep(INIT_SIZE_TYPE, long InitK) 
527   { k = MaxK = -1; tbl = 0; NumPrimes = 0; SetSize(InitK); }
528   ~FFTRep();
529};
530
531
532void ToFFTRep(FFTRep& y, const ZZ_pX& x, long k, long lo, long hi);
533// computes an n = 2^k point convolution of x[lo..hi].
534
535inline void ToFFTRep(FFTRep& y, const ZZ_pX& x, long k)
536
537   { ToFFTRep(y, x, k, 0, deg(x)); }
538
539void RevToFFTRep(FFTRep& y, const vec_ZZ_p& x,
540                 long k, long lo, long hi, long offset);
541// computes an n = 2^k point convolution of X^offset*x[lo..hi] mod X^n-1
542// using "inverted" evaluation points.
543
544
545void FromFFTRep(ZZ_pX& x, FFTRep& y, long lo, long hi);
546// converts from FFT-representation to coefficient representation
547// only the coefficients lo..hi are computed
548// NOTE: this version destroys the data in y
549
550// non-destructive versions of the above
551
552void NDFromFFTRep(ZZ_pX& x, const FFTRep& y, long lo, long hi, FFTRep& temp);
553void NDFromFFTRep(ZZ_pX& x, const FFTRep& y, long lo, long hi);
554
555void RevFromFFTRep(vec_ZZ_p& x, FFTRep& y, long lo, long hi);
556
557   // converts from FFT-representation to coefficient representation
558   // using "inverted" evaluation points.
559   // only the coefficients lo..hi are computed
560
561
562
563
564void FromFFTRep(ZZ_p* x, FFTRep& y, long lo, long hi);
565// convert out coefficients lo..hi of y, store result in x.
566// no normalization is done.
567
568
569// direct manipulation of FFT reps
570
571void mul(FFTRep& z, const FFTRep& x, const FFTRep& y);
572void sub(FFTRep& z, const FFTRep& x, const FFTRep& y);
573void add(FFTRep& z, const FFTRep& x, const FFTRep& y);
574
575void reduce(FFTRep& x, const FFTRep& a, long k);
576// reduces a 2^l point FFT-rep to a 2^k point FFT-rep
577
578void AddExpand(FFTRep& x, const FFTRep& a);
579//  x = x + (an "expanded" version of a)
580
581
582
583
584
585// This data structure holds unconvoluted modular representations
586// of polynomials
587
588class ZZ_pXModRep {
589private:
590   ZZ_pXModRep(const ZZ_pXModRep&); // disabled
591   void operator=(const ZZ_pXModRep&); // disabled
592
593public:
594   long n;
595   long MaxN;
596   long **tbl;
597   long NumPrimes; 
598
599   void SetSize(long NewN);
600
601   ZZ_pXModRep() { n = MaxN = 0; tbl = 0; NumPrimes = 0; }
602   ZZ_pXModRep(INIT_SIZE_TYPE, long k) 
603   { n = MaxN = 0; tbl = 0; NumPrimes = 0; SetSize(k); }
604   ~ZZ_pXModRep();
605};
606
607
608void ToZZ_pXModRep(ZZ_pXModRep& x, const ZZ_pX& a, long lo, long hi);
609
610void ToFFTRep(FFTRep& x, const ZZ_pXModRep& a, long k, long lo, long hi);
611// converts coefficients lo..hi to a 2^k-point FFTRep.
612// must have hi-lo+1 < 2^k
613
614
615
616
617
618/*************************************************************
619
620                      Division
621
622**************************************************************/
623
624void DivRem(ZZ_pX& q, ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b);
625// q = a/b, r = a%b
626
627void div(ZZ_pX& q, const ZZ_pX& a, const ZZ_pX& b);
628// q = a/b
629
630void div(ZZ_pX& q, const ZZ_pX& a, const ZZ_p& b);
631void div(ZZ_pX& q, const ZZ_pX& a, long b);
632
633
634void rem(ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b);
635// r = a%b
636
637long divide(ZZ_pX& q, const ZZ_pX& a, const ZZ_pX& b);
638// if b | a, sets q = a/b and returns 1; otherwise returns 0
639
640long divide(const ZZ_pX& a, const ZZ_pX& b);
641// if b | a, sets q = a/b and returns 1; otherwise returns 0
642
643void InvTrunc(ZZ_pX& x, const ZZ_pX& a, long m);
644// computes x = a^{-1} % X^m
645// constant term must be non-zero
646
647inline ZZ_pX InvTrunc(const ZZ_pX& a, long m)
648   { ZZ_pX x; InvTrunc(x, a, m); NTL_OPT_RETURN(ZZ_pX, x); }
649
650
651
652// These always use "classical" arithmetic
653void PlainDivRem(ZZ_pX& q, ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b);
654void PlainDiv(ZZ_pX& q, const ZZ_pX& a, const ZZ_pX& b);
655void PlainRem(ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b);
656
657void PlainRem(ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b, ZZVec& tmp);
658void PlainDivRem(ZZ_pX& q, ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b,
659                 ZZVec& tmp);
660
661
662// These always use FFT arithmetic
663void FFTDivRem(ZZ_pX& q, ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b);
664void FFTDiv(ZZ_pX& q, const ZZ_pX& a, const ZZ_pX& b);
665void FFTRem(ZZ_pX& r, const ZZ_pX& a, const ZZ_pX& b);
666
667void PlainInvTrunc(ZZ_pX& x, const ZZ_pX& a, long m);
668// always uses "classical" algorithm
669// ALIAS RESTRICTION: input may not alias output
670
671void NewtonInvTrunc(ZZ_pX& x, const ZZ_pX& a, long m);
672// uses a Newton Iteration with the FFT.
673// ALIAS RESTRICTION: input may not alias output
674
675
676inline ZZ_pX operator/(const ZZ_pX& a, const ZZ_pX& b)
677   { ZZ_pX x; div(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
678
679inline ZZ_pX operator/(const ZZ_pX& a, const ZZ_p& b)
680   { ZZ_pX x; div(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
681
682inline ZZ_pX operator/(const ZZ_pX& a, long b)
683   { ZZ_pX x; div(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
684
685inline ZZ_pX& operator/=(ZZ_pX& x, const ZZ_p& b)
686   { div(x, x, b); return x; }
687
688inline ZZ_pX& operator/=(ZZ_pX& x, long b)
689   { div(x, x, b); return x; }
690
691inline ZZ_pX& operator/=(ZZ_pX& x, const ZZ_pX& b)
692   { div(x, x, b); return x; }
693
694
695inline ZZ_pX operator%(const ZZ_pX& a, const ZZ_pX& b)
696   { ZZ_pX x; rem(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
697
698inline ZZ_pX& operator%=(ZZ_pX& x, const ZZ_pX& b)
699   { rem(x, x, b); return x; }
700
701
702/***********************************************************
703
704                         GCD's
705
706************************************************************/
707
708
709void GCD(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
710// x = GCD(a, b),  x is always monic (or zero if a==b==0).
711
712inline ZZ_pX GCD(const ZZ_pX& a, const ZZ_pX& b)
713   { ZZ_pX x; GCD(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
714
715void XGCD(ZZ_pX& d, ZZ_pX& s, ZZ_pX& t, const ZZ_pX& a, const ZZ_pX& b);
716// d = gcd(a,b), a s + b t = d
717
718void PlainXGCD(ZZ_pX& d, ZZ_pX& s, ZZ_pX& t, const ZZ_pX& a, const ZZ_pX& b);
719// same as above, but uses classical algorithm
720
721
722void PlainGCD(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b);
723// always uses "cdlassical" arithmetic
724
725
726class ZZ_pXMatrix {
727private:
728
729   ZZ_pXMatrix(const ZZ_pXMatrix&);  // disable
730   ZZ_pX elts[2][2];
731
732public:
733
734   ZZ_pXMatrix() { }
735   ~ZZ_pXMatrix() { }
736
737   void operator=(const ZZ_pXMatrix&);
738   ZZ_pX& operator() (long i, long j) { return elts[i][j]; }
739   const ZZ_pX& operator() (long i, long j) const { return elts[i][j]; }
740};
741
742
743void HalfGCD(ZZ_pXMatrix& M_out, const ZZ_pX& U, const ZZ_pX& V, long d_red);
744// deg(U) > deg(V),   1 <= d_red <= deg(U)+1.
745//
746// This computes a 2 x 2 polynomial matrix M_out such that
747//    M_out * (U, V)^T = (U', V')^T,
748// where U', V' are consecutive polynomials in the Euclidean remainder
749// sequence of U, V, and V' is the polynomial of highest degree
750// satisfying deg(V') <= deg(U) - d_red.
751
752void XHalfGCD(ZZ_pXMatrix& M_out, ZZ_pX& U, ZZ_pX& V, long d_red);
753
754// same as above, except that U is replaced by U', and V by V'
755
756
757/*************************************************************
758
759             Modular Arithmetic without pre-conditioning
760
761**************************************************************/
762
763// arithmetic mod f.
764// all inputs and outputs are polynomials of degree less than deg(f).
765// ASSUMPTION: f is assumed monic, and deg(f) > 0.
766// NOTE: if you want to do many computations with a fixed f,
767//       use the ZZ_pXModulus data structure and associated routines below.
768
769
770
771void MulMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b, const ZZ_pX& f);
772// x = (a * b) % f
773
774inline ZZ_pX MulMod(const ZZ_pX& a, const ZZ_pX& b, const ZZ_pX& f)
775   { ZZ_pX x; MulMod(x, a, b, f); NTL_OPT_RETURN(ZZ_pX, x); }
776
777void SqrMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& f);
778// x = a^2 % f
779
780inline ZZ_pX SqrMod(const ZZ_pX& a,  const ZZ_pX& f)
781   { ZZ_pX x; SqrMod(x, a, f); NTL_OPT_RETURN(ZZ_pX, x); }
782
783void MulByXMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& f);
784// x = (a * X) mod f
785
786inline ZZ_pX MulByXMod(const ZZ_pX& a,  const ZZ_pX& f)
787   { ZZ_pX x; MulByXMod(x, a, f); NTL_OPT_RETURN(ZZ_pX, x); }
788
789
790
791void InvMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& f);
792// x = a^{-1} % f, error is a is not invertible
793
794inline ZZ_pX InvMod(const ZZ_pX& a,  const ZZ_pX& f)
795   { ZZ_pX x; InvMod(x, a, f); NTL_OPT_RETURN(ZZ_pX, x); }
796
797long InvModStatus(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& f);
798// if (a, f) = 1, returns 0 and sets x = a^{-1} % f
799// otherwise, returns 1 and sets x = (a, f)
800
801
802
803
804/******************************************************************
805
806        Modular Arithmetic with Pre-conditioning
807
808*******************************************************************/
809
810
811// If you need to do a lot of arithmetic modulo a fixed f,
812// build ZZ_pXModulus F for f.  This pre-computes information about f
813// that speeds up the computation a great deal.
814
815
816class ZZ_pXModulus {
817public:
818   ZZ_pXModulus() : UseFFT(0), n(-1)  { }
819   ~ZZ_pXModulus() { }
820   
821
822   // the following members may become private in future
823   ZZ_pX f;   // the modulus
824   long UseFFT;// flag indicating whether FFT should be used.
825   long n;     // n = deg(f)
826   long k;     // least k s/t 2^k >= n
827   long l;     // least l s/t 2^l >= 2n-3
828   FFTRep FRep; // 2^k point rep of f
829                // H = rev((rev(f))^{-1} rem X^{n-1})
830   FFTRep HRep; // 2^l point rep of H
831   vec_ZZ_p tracevec;  // mutable
832
833   // but these will remain public
834   ZZ_pXModulus(const ZZ_pX& ff);
835
836   const ZZ_pX& val() const { return f; }
837   operator const ZZ_pX& () const { return f; }
838
839};
840
841inline long deg(const ZZ_pXModulus& F) { return F.n; }
842
843void build(ZZ_pXModulus& F, const ZZ_pX& f);
844// deg(f) > 0.
845
846
847void rem21(ZZ_pX& x, const ZZ_pX& a, const ZZ_pXModulus& F);
848// x = a % f
849// deg(a) <= 2(n-1), where n = F.n = deg(f)
850
851void rem(ZZ_pX& x, const ZZ_pX& a, const ZZ_pXModulus& F);
852// x = a % f, no restrictions on deg(a);  makes repeated calls to rem21
853
854inline ZZ_pX operator%(const ZZ_pX& a, const ZZ_pXModulus& F)
855   { ZZ_pX x; rem(x, a, F); NTL_OPT_RETURN(ZZ_pX, x); }
856
857inline ZZ_pX& operator%=(ZZ_pX& x, const ZZ_pXModulus& F)
858   { rem(x, x, F); return x; } 
859
860void DivRem(ZZ_pX& q, ZZ_pX& r, const ZZ_pX& a, const ZZ_pXModulus& F);
861
862void div(ZZ_pX& q, const ZZ_pX& a, const ZZ_pXModulus& F);
863
864inline ZZ_pX operator/(const ZZ_pX& a, const ZZ_pXModulus& F)
865   { ZZ_pX x; div(x, a, F); NTL_OPT_RETURN(ZZ_pX, x); }
866
867inline ZZ_pX& operator/=(ZZ_pX& x, const ZZ_pXModulus& F)
868   { div(x, x, F); return x; } 
869
870void MulMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pX& b, const ZZ_pXModulus& F);
871// x = (a * b) % f
872// deg(a), deg(b) < n
873
874inline ZZ_pX MulMod(const ZZ_pX& a, const ZZ_pX& b, const ZZ_pXModulus& F)
875   { ZZ_pX x; MulMod(x, a, b, F); NTL_OPT_RETURN(ZZ_pX, x); }
876
877void SqrMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pXModulus& F);
878// x = a^2 % f
879// deg(a) < n
880
881inline ZZ_pX SqrMod(const ZZ_pX& a, const ZZ_pXModulus& F)
882   { ZZ_pX x; SqrMod(x, a, F); NTL_OPT_RETURN(ZZ_pX, x); }
883
884void PowerMod(ZZ_pX& x, const ZZ_pX& a, const ZZ& e, const ZZ_pXModulus& F);
885// x = a^e % f, e >= 0
886
887inline ZZ_pX PowerMod(const ZZ_pX& a, const ZZ& e, const ZZ_pXModulus& F)
888   { ZZ_pX x; PowerMod(x, a, e, F); NTL_OPT_RETURN(ZZ_pX, x); }
889
890inline void PowerMod(ZZ_pX& x, const ZZ_pX& a, long e, const ZZ_pXModulus& F)
891   { PowerMod(x, a, ZZ_expo(e), F); }
892
893inline ZZ_pX PowerMod(const ZZ_pX& a, long e, const ZZ_pXModulus& F)
894   { ZZ_pX x; PowerMod(x, a, e, F); NTL_OPT_RETURN(ZZ_pX, x); }
895
896
897
898void PowerXMod(ZZ_pX& x, const ZZ& e, const ZZ_pXModulus& F);
899// x = X^e % f, e >= 0
900
901inline ZZ_pX PowerXMod(const ZZ& e, const ZZ_pXModulus& F)
902   { ZZ_pX x; PowerXMod(x, e, F); NTL_OPT_RETURN(ZZ_pX, x); }
903
904inline void PowerXMod(ZZ_pX& x, long e, const ZZ_pXModulus& F)
905   { PowerXMod(x, ZZ_expo(e), F); }
906
907inline ZZ_pX PowerXMod(long e, const ZZ_pXModulus& F)
908   { ZZ_pX x; PowerXMod(x, e, F); NTL_OPT_RETURN(ZZ_pX, x); }
909
910void PowerXPlusAMod(ZZ_pX& x, const ZZ_p& a, const ZZ& e, const ZZ_pXModulus& F);
911// x = (X + a)^e % f, e >= 0
912
913inline ZZ_pX PowerXPlusAMod(const ZZ_p& a, const ZZ& e, const ZZ_pXModulus& F)
914   { ZZ_pX x; PowerXPlusAMod(x, a, e, F); NTL_OPT_RETURN(ZZ_pX, x); }
915
916inline void PowerXPlusAMod(ZZ_pX& x, const ZZ_p& a, long e, const ZZ_pXModulus& F)
917   { PowerXPlusAMod(x, a, ZZ_expo(e), F); }
918
919
920inline ZZ_pX PowerXPlusAMod(const ZZ_p& a, long e, const ZZ_pXModulus& F)
921   { ZZ_pX x; PowerXPlusAMod(x, a, e, F); NTL_OPT_RETURN(ZZ_pX, x); }
922
923// If you need to compute a * b % f for a fixed b, but for many a's
924// (for example, computing powers of b modulo f), it is
925// much more efficient to first build a ZZ_pXMultiplier B for b,
926// and then use the routine below.
927
928class ZZ_pXMultiplier {
929public:
930   ZZ_pXMultiplier() : UseFFT(0)  { }
931   ZZ_pXMultiplier(const ZZ_pX& b, const ZZ_pXModulus& F);
932
933   ~ZZ_pXMultiplier() { }
934
935
936   // the following members may become private in the future
937   ZZ_pX b;   
938   long UseFFT;
939   FFTRep B1; 
940   FFTRep B2; 
941
942   // but this will remain public
943   const ZZ_pX& val() const { return b; }
944
945};
946
947void build(ZZ_pXMultiplier& B, const ZZ_pX& b, const ZZ_pXModulus& F);
948
949void MulMod(ZZ_pX& x, const ZZ_pX& a, const ZZ_pXMultiplier& B,
950                                      const ZZ_pXModulus& F);
951
952inline ZZ_pX MulMod(const ZZ_pX& a, const ZZ_pXMultiplier& B,
953                                          const ZZ_pXModulus& F)
954   { ZZ_pX x; MulMod(x, a, B, F); NTL_OPT_RETURN(ZZ_pX, x); }
955
956// x = (a * b) % f
957
958
959/*******************************************************
960
961              Evaluation and related problems
962
963********************************************************/
964
965
966void BuildFromRoots(ZZ_pX& x, const vec_ZZ_p& a);
967// computes the polynomial (X-a[0]) ... (X-a[n-1]), where n = a.length()
968
969inline ZZ_pX BuildFromRoots(const vec_ZZ_p& a)
970   { ZZ_pX x; BuildFromRoots(x, a); NTL_OPT_RETURN(ZZ_pX, x); }
971
972
973void eval(ZZ_p& b, const ZZ_pX& f, const ZZ_p& a);
974// b = f(a)
975
976inline ZZ_p eval(const ZZ_pX& f, const ZZ_p& a)
977   { ZZ_p x; eval(x, f, a); NTL_OPT_RETURN(ZZ_p, x); }
978
979void eval(vec_ZZ_p& b, const ZZ_pX& f, const vec_ZZ_p& a);
980//  b[i] = f(a[i])
981
982inline vec_ZZ_p eval(const ZZ_pX& f, const vec_ZZ_p& a)
983   { vec_ZZ_p x; eval(x, f, a); NTL_OPT_RETURN(vec_ZZ_p, x); }
984
985
986void interpolate(ZZ_pX& f, const vec_ZZ_p& a, const vec_ZZ_p& b);
987// computes f such that f(a[i]) = b[i]
988
989inline ZZ_pX interpolate(const vec_ZZ_p& a, const vec_ZZ_p& b)
990   { ZZ_pX x; interpolate(x, a, b); NTL_OPT_RETURN(ZZ_pX, x); }
991
992
993/*****************************************************************
994
995                       vectors of ZZ_pX's
996
997*****************************************************************/
998
999NTL_vector_decl(ZZ_pX,vec_ZZ_pX)
1000
1001NTL_eq_vector_decl(ZZ_pX,vec_ZZ_pX)
1002
1003NTL_io_vector_decl(ZZ_pX,vec_ZZ_pX)
1004
1005
1006
1007/**********************************************************
1008
1009         Modular Composition and Minimal Polynomials
1010
1011***********************************************************/
1012
1013
1014// algorithms for computing g(h) mod f
1015
1016
1017
1018void CompMod(ZZ_pX& x, const ZZ_pX& g, const ZZ_pX& h, const ZZ_pXModulus& F);
1019// x = g(h) mod f
1020
1021inline ZZ_pX CompMod(const ZZ_pX& g, const ZZ_pX& h, 
1022                           const ZZ_pXModulus& F)
1023   { ZZ_pX x; CompMod(x, g, h, F); NTL_OPT_RETURN(ZZ_pX, x); }
1024
1025void Comp2Mod(ZZ_pX& x1, ZZ_pX& x2, const ZZ_pX& g1, const ZZ_pX& g2,
1026              const ZZ_pX& h, const ZZ_pXModulus& F);
1027// xi = gi(h) mod f (i=1,2)
1028
1029void Comp3Mod(ZZ_pX& x1, ZZ_pX& x2, ZZ_pX& x3, 
1030              const ZZ_pX& g1, const ZZ_pX& g2, const ZZ_pX& g3,
1031              const ZZ_pX& h, const ZZ_pXModulus& F);
1032// xi = gi(h) mod f (i=1..3)
1033
1034
1035
1036// The routine build (see below) which is implicitly called
1037// by the various compose and UpdateMap routines builds a table
1038// of polynomials. 
1039// If ZZ_pXArgBound > 0, then the table is limited in
1040// size to approximamtely that many KB.
1041// If ZZ_pXArgBound <= 0, then it is ignored, and space is allocated
1042// so as to maximize speed.
1043// Initially, ZZ_pXArgBound = 0.
1044
1045
1046// If a single h is going to be used with many g's
1047// then you should build a ZZ_pXArgument for h,
1048// and then use the compose routine below.
1049// build computes and stores h, h^2, ..., h^m mod f.
1050// After this pre-computation, composing a polynomial of degree
1051// roughly n with h takes n/m multiplies mod f, plus n^2
1052// scalar multiplies.
1053// Thus, increasing m increases the space requirement and the pre-computation
1054// time, but reduces the composition time.
1055// If ZZ_pXArgBound > 0, a table of size less than m may be built.
1056
1057struct ZZ_pXArgument {
1058   vec_ZZ_pX H;
1059};
1060
1061extern long ZZ_pXArgBound;
1062
1063
1064void build(ZZ_pXArgument& H, const ZZ_pX& h, const ZZ_pXModulus& F, long m);
1065
1066// m must be > 0, otherwise an error is raised
1067
1068void CompMod(ZZ_pX& x, const ZZ_pX& g, const ZZ_pXArgument& H, 
1069             const ZZ_pXModulus& F);
1070
1071inline ZZ_pX
1072CompMod(const ZZ_pX& g, const ZZ_pXArgument& H, const ZZ_pXModulus& F)
1073   { ZZ_pX x; CompMod(x, g, H, F); NTL_OPT_RETURN(ZZ_pX, x); }
1074
1075
1076#ifndef NTL_TRANSITION
1077
1078void UpdateMap(vec_ZZ_p& x, const vec_ZZ_p& a,
1079               const ZZ_pXMultiplier& B, const ZZ_pXModulus& F);
1080
1081inline vec_ZZ_p
1082UpdateMap(const vec_ZZ_p& a, 
1083          const ZZ_pXMultiplier& B, const ZZ_pXModulus& F)
1084   { vec_ZZ_p x; UpdateMap(x, a, B, F); 
1085     NTL_OPT_RETURN(vec_ZZ_p, x); }
1086
1087#endif
1088
1089
1090/* computes (a, b), (a, (b*X)%f), ..., (a, (b*X^{n-1})%f),
1091   where ( , ) denotes the vector inner product.
1092
1093   This is really a "transposed" MulMod by B.
1094*/
1095
1096void PlainUpdateMap(vec_ZZ_p& x, const vec_ZZ_p& a,
1097                    const ZZ_pX& b, const ZZ_pX& f);
1098
1099// same as above, but uses only classical arithmetic
1100
1101
1102void ProjectPowers(vec_ZZ_p& x, const vec_ZZ_p& a, long k,
1103                   const ZZ_pX& h, const ZZ_pXModulus& F);
1104
1105inline vec_ZZ_p ProjectPowers(const vec_ZZ_p& a, long k,
1106                   const ZZ_pX& h, const ZZ_pXModulus& F)
1107{
1108   vec_ZZ_p x;
1109   ProjectPowers(x, a, k, h, F);
1110   NTL_OPT_RETURN(vec_ZZ_p, x);
1111}
1112
1113
1114// computes (a, 1), (a, h), ..., (a, h^{k-1} % f)
1115// this is really a "transposed" compose.
1116
1117void ProjectPowers(vec_ZZ_p& x, const vec_ZZ_p& a, long k,
1118                   const ZZ_pXArgument& H, const ZZ_pXModulus& F);
1119
1120inline vec_ZZ_p ProjectPowers(const vec_ZZ_p& a, long k,
1121                   const ZZ_pXArgument& H, const ZZ_pXModulus& F)
1122{
1123   vec_ZZ_p x;
1124   ProjectPowers(x, a, k, H, F);
1125   NTL_OPT_RETURN(vec_ZZ_p, x);
1126}
1127
1128// same as above, but uses a pre-computed ZZ_pXArgument
1129
1130inline void project(ZZ_p& x, const vec_ZZ_p& a, const ZZ_pX& b)
1131   { InnerProduct(x, a, b.rep); }
1132
1133inline ZZ_p project(const vec_ZZ_p& a, const ZZ_pX& b)
1134   {  ZZ_p x; project(x, a, b); NTL_OPT_RETURN(ZZ_p, x); }
1135
1136void MinPolySeq(ZZ_pX& h, const vec_ZZ_p& a, long m);
1137// computes the minimum polynomial of a linealy generated sequence;
1138// m is a bound on the degree of the polynomial;
1139// required: a.length() >= 2*m
1140
1141inline ZZ_pX MinPolySeq(const vec_ZZ_p& a, long m)
1142   { ZZ_pX x; MinPolySeq(x, a, m); NTL_OPT_RETURN(ZZ_pX, x); }
1143
1144void ProbMinPolyMod(ZZ_pX& h, const ZZ_pX& g, const ZZ_pXModulus& F, long m);
1145
1146inline ZZ_pX ProbMinPolyMod(const ZZ_pX& g, const ZZ_pXModulus& F, long m)
1147   { ZZ_pX x; ProbMinPolyMod(x, g, F, m); NTL_OPT_RETURN(ZZ_pX, x); }
1148
1149
1150inline void ProbMinPolyMod(ZZ_pX& h, const ZZ_pX& g, const ZZ_pXModulus& F)
1151   { ProbMinPolyMod(h, g, F, F.n); }
1152
1153inline ZZ_pX ProbMinPolyMod(const ZZ_pX& g, const ZZ_pXModulus& F)
1154   { ZZ_pX x; ProbMinPolyMod(x, g, F); NTL_OPT_RETURN(ZZ_pX, x); }
1155
1156
1157// computes the monic minimal polynomial if (g mod f).
1158// m = a bound on the degree of the minimal polynomial.
1159// If this argument is not supplied, it defaults to deg(f).
1160// The algorithm is probabilistic, always returns a divisor of
1161// the minimal polynomial, and returns a proper divisor with
1162// probability at most m/p.
1163
1164void MinPolyMod(ZZ_pX& h, const ZZ_pX& g, const ZZ_pXModulus& F, long m);
1165
1166inline ZZ_pX MinPolyMod(const ZZ_pX& g, const ZZ_pXModulus& F, long m)
1167   { ZZ_pX x; MinPolyMod(x, g, F, m); NTL_OPT_RETURN(ZZ_pX, x); }
1168
1169inline void MinPolyMod(ZZ_pX& h, const ZZ_pX& g, const ZZ_pXModulus& F)
1170   { MinPolyMod(h, g, F, F.n); }
1171
1172inline ZZ_pX MinPolyMod(const ZZ_pX& g, const ZZ_pXModulus& F)
1173   { ZZ_pX x; MinPolyMod(x, g, F); NTL_OPT_RETURN(ZZ_pX, x); }
1174
1175// same as above, but guarantees that result is correct
1176
1177void IrredPolyMod(ZZ_pX& h, const ZZ_pX& g, const ZZ_pXModulus& F, long m);
1178
1179inline ZZ_pX IrredPolyMod(const ZZ_pX& g, const ZZ_pXModulus& F, long m)
1180   { ZZ_pX x; IrredPolyMod(x, g, F, m); NTL_OPT_RETURN(ZZ_pX, x); }
1181
1182
1183inline void IrredPolyMod(ZZ_pX& h, const ZZ_pX& g, const ZZ_pXModulus& F)
1184   { IrredPolyMod(h, g, F, F.n); }
1185
1186inline ZZ_pX IrredPolyMod(const ZZ_pX& g, const ZZ_pXModulus& F)
1187   { ZZ_pX x; IrredPolyMod(x, g, F); NTL_OPT_RETURN(ZZ_pX, x); }
1188
1189// same as above, but assumes that f is irreducible,
1190// or at least that the minimal poly of g is itself irreducible.
1191// The algorithm is deterministic (and is always correct).
1192
1193/*****************************************************************
1194
1195                   Traces, norms, resultants
1196
1197******************************************************************/
1198
1199void TraceVec(vec_ZZ_p& S, const ZZ_pX& f);
1200
1201inline vec_ZZ_p TraceVec(const ZZ_pX& f)
1202   { vec_ZZ_p x; TraceVec(x, f); NTL_OPT_RETURN(vec_ZZ_p, x); }
1203
1204void FastTraceVec(vec_ZZ_p& S, const ZZ_pX& f);
1205void PlainTraceVec(vec_ZZ_p& S, const ZZ_pX& f);
1206
1207void TraceMod(ZZ_p& x, const ZZ_pX& a, const ZZ_pXModulus& F);
1208
1209inline ZZ_p TraceMod(const ZZ_pX& a, const ZZ_pXModulus& F)
1210   { ZZ_p x; TraceMod(x, a, F); NTL_OPT_RETURN(ZZ_p, x); }
1211
1212void TraceMod(ZZ_p& x, const ZZ_pX& a, const ZZ_pX& f);
1213
1214inline ZZ_p TraceMod(const ZZ_pX& a, const ZZ_pX& f)
1215   { ZZ_p x; TraceMod(x, a, f); NTL_OPT_RETURN(ZZ_p, x); }
1216
1217
1218
1219void ComputeTraceVec(const ZZ_pXModulus& F);
1220
1221
1222void NormMod(ZZ_p& x, const ZZ_pX& a, const ZZ_pX& f);
1223
1224inline ZZ_p NormMod(const ZZ_pX& a, const ZZ_pX& f)
1225   { ZZ_p x; NormMod(x, a, f); NTL_OPT_RETURN(ZZ_p, x); }
1226
1227void resultant(ZZ_p& rres, const ZZ_pX& a, const ZZ_pX& b);
1228
1229inline ZZ_p resultant(const ZZ_pX& a, const ZZ_pX& b)
1230   { ZZ_p x; resultant(x, a, b); NTL_OPT_RETURN(ZZ_p, x); }
1231
1232void CharPolyMod(ZZ_pX& g, const ZZ_pX& a, const ZZ_pX& f);
1233// g = char poly of (a mod f)
1234
1235inline ZZ_pX CharPolyMod(const ZZ_pX& a, const ZZ_pX& f)
1236   { ZZ_pX x; CharPolyMod(x, a, f); NTL_OPT_RETURN(ZZ_pX, x); }
1237
1238
1239
1240NTL_CLOSE_NNS
1241
1242#endif
Note: See TracBrowser for help on using the repository browser.