-
Notifications
You must be signed in to change notification settings - Fork 1
/
quads.cc
1173 lines (1064 loc) · 35 KB
/
quads.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// FILE QUADS.CC
#define testbezout // define this to turn on self-verification of bezout
#include <iostream>
#include "intprocs.h"
#include "quads.h"
#include "primes.h"
#include "geometry.h"
#include "swan.h"
//Declare static data members of class Quad (except SD which is declared in geometry.cc)
long Quad::d;
INT Quad::disc;
INT Quad::absdisc;
vector<INT> Quad::prime_disc_factors;
vector<INT> Quad::all_disc_factors;
long Quad::t;
INT Quad::n;
char Quad::name;
long Quad::maxnorm;
int Quad::nunits;
int Quad::is_Euclidean;
int Quad::class_number;
int Quad::class_group_2_rank;
Quad Quad::w;
Quad Quad::zero;
Quad Quad::one;
vector<Quad> Quad::shifts_by_one;
int Quad::geometry_initialised;
//Primes
vector<Quad> quadprimes; //Initialised by initquadprimes, see below
long nquadprimes; //The number of them.
vector<Quad> quadunits, squareunits;
Quad fundunit;
vector<long> euclidean_fields = {1,2,3,7,11};
// fields for which geometry is defined (all 305 with |disc|<1000 *except* 246 (984) and 249 (996), which include all h<=3)
vector<long> valid_fields = {1, 2, 3, 5, 6, 7,
10, 11, 13, 14, 15, 17, 19, 21, 22, 23, 26, 29, 30, 31, 33, 34, 35, 37, 38, 39,
41, 42, 43, 46, 47, 51, 53, 55, 57, 58, 59, 61, 62, 65, 66, 67, 69, 70, 71, 73,
74, 77, 78, 79, 82, 83, 85, 86, 87, 89, 91, 93, 94, 95, 97,
101, 102, 103, 105, 106, 107, 109, 110, 111, 113, 114, 115, 118, 119, 122, 123,
127, 129, 130, 131, 133, 134, 137, 138, 139, 141, 142, 143, 145, 146, 149, 151,
154, 155, 157, 158, 159, 161, 163, 165, 166, 167, 170, 173, 174, 177, 178, 179,
181, 182, 183, 185, 186, 187, 190, 191, 193, 194, 195, 197, 199,
201, 202, 203, 205, 206, 209, 210, 211, 213, 214, 215, 217, 218, 219, 221, 222,
223, 226, 227, 229, 230, 231, 233, 235, 237, 238, 239, 241, 246, 247, 249, 251,
253, 254, 255, 257, 258, 259, 262, 263, 265, 266, 269, 267, 271, 273, 274, 277,
278, 281, 282, 283, 285, 286, 287, 290, 291, 293, 295, 298, 299,
301, 302, 303, 305, 307, 309, 310, 311, 314, 313, 317, 318, 319, 322, 321, 323,
326, 327, 329, 330, 331, 334, 335, 337, 339, 341, 345, 346, 347, 349, 353, 354,
355, 357, 358, 359, 362, 365, 366, 367, 370, 371, 373, 374, 377, 379, 381, 382,
383, 385, 386, 389, 390, 391, 393, 394, 395, 397, 398, 399,
401, 402, 403, 406, 407, 409, 410, 411, 413, 415, 417, 418, 419, 421, 422, 426,
427, 429, 430, 431, 433, 434, 435, 437, 438, 439, 442, 443, 445, 446, 447, 449,
451, 453, 454, 455, 457, 458, 461, 462, 463, 465, 466, 467, 469, 470, 471, 473,
474, 478, 479, 481, 482, 483, 485, 487, 489, 491, 493, 494, 497, 498, 499,
501, 502, 503, 505, 506, 509, 510, 511, 514, 515, 517, 518, 519, 521, 523, 526, 527,
530, 535, 543, 547, 551, 553, 555, 559, 562, 563, 571, 579, 583, 587, 591, 595, 599,
607, 611, 615, 619, 623, 627, 631, 635, 643, 647, 651, 655, 659, 663, 667, 671,
679, 683, 687, 691, 695, 699,
703, 707, 715, 719, 723, 727, 731, 739, 743, 751, 755, 759, 763, 767, 771, 779,
787, 791, 795, 799,
803, 807, 811, 815, 823, 827, 831, 835, 839, 843, 851, 859, 863, 871, 879, 883,
887, 895, 899,
903, 907, 911, 915, 919, 923, 935, 939, 943, 947, 951, 955, 959, 967, 971, 979,
983, 987, 991, 995,
1003, 1007, 1011, 1015, 1019, 1023, 1027, 1031, 1039, 1043, 1047, 1051, 1055, 1059,
1063, 1067, 1079, 1087, 1091, 1095, 1099,
1103, 1111, 1115, 1119, 1123, 1131, 1135, 1139, 1147, 1151, 1155, 1159, 1163, 1167,
1171, 1187, 1191, 1195, 1199,
1203, 1207, 1211, 1219, 1223, 1227, 1231, 1235, 1239, 1243, 1247, 1255, 1259, 1263,
1267, 1271, 1279, 1283, 1291, 1295, 1299,
1303, 1307, 1311, 1315, 1319, 1327, 1335, 1339, 1343, 1347, 1351, 1355, 1363, 1367,
1371, 1379, 1383, 1387, 1391, 1399,
1403, 1407, 1411, 1415, 1419, 1423, 1427, 1435, 1439, 1443, 1447, 1451, 1455, 1459,
1463, 1471, 1479, 1483, 1487, 1491, 1495, 1499,
1507, 1511, 1515, 1523, 1527, 1531, 1535, 1543, 1547, 1551, 1555, 1559, 1563, 1567,
1571, 1579, 1583, 1591, 1595, 1599,
1603, 1607, 1615, 1619, 1623, 1627, 1631, 1635, 1639, 1643, 1651, 1655, 1659, 1663,
1667, 1671, 1679, 1687, 1691, 1695, 1699,
1703, 1707, 1711, 1723, 1727, 1731, 1735, 1739, 1743, 1751, 1747, 1759, 1763, 1767,
1771, 1779, 1783, 1787, 1795, 1799,
1803, 1807, 1819, 1811, 1823, 1831, 1835, 1839, 1843, 1847, 1851, 1855, 1867, 1871,
1879, 1883, 1887, 1891, 1895,
1903, 1907, 1915, 1919, 1923, 1927, 1931, 1939, 1943, 1947, 1951, 1955, 1959, 1963,
1967, 1979, 1983, 1987, 1991, 1995, 1999,
2003, 2011, 2015, 2019, 2027, 2031, 2035, 2039, 2047, 2051, 2055, 2059, 2063, 2067,
2071, 2083, 2087, 2091, 2095, 2099,
2103, 2111, 2119, 2123, 2127, 2131, 2135, 2139, 2143, 2147, 2155, 2159, 2167, 2171,
2179, 2183, 2191, 2195,
2203, 2207, 2211, 2215, 2219, 2227, 2231, 2235, 2239, 2243, 2247, 2251, 2255, 2263,
2267, 2271, 2279, 2283, 2287, 2291,
2347,
2683};
vector<long> class_number_one_fields = {1,2,3,7,11,19,43,67,163}; // 9
vector<long> class_number_two_fields = {5,6,10,13,15,22,35,37,51,58,91,115,123,187,235,267,403,427}; // 18
vector<long> class_number_three_fields = {23,31,59,83,107,139,211,283,307,331,379,499,547,643,883,907}; // 16
vector<long> class_number_four_fields = {14,17,21,30,33,34,39,42,46,55,57,70,73,78,82,85,93,97,102,130,
133,142,155,177,190,193,195,203,219,253,259,291,323,355,435,483,
555,595,627,667,715,723,763,795,955,1003,1027,1227,1243,1387,1411,
1435,1507, 1555}; // 54
vector<long> class_number_five_fields = {47,79,103,127,131,179,227,347,443,523,571,619,683,691,739,787,947,
1051,1123,1723,1747,1867,2203,2347, 2683}; // 25
vector<long> valid_field_discs(long max_disc)
{
vector<long> discs;
for(auto d : valid_fields)
{
long D = (d%4==3? d : 4*d);
if ((max_disc==0) || (D<=max_disc))
discs.push_back(D);
}
::sort(discs.begin(), discs.end());
return discs;
}
int check_field(long d, vector<long> fields)
{
return (std::find(fields.begin(), fields.end(), d) != fields.end());
}
// declaration of "extern" functions declared in quads.h:
Quad (*mult)(const Quad& a, const Quad& b);
Quad (*qdivi)(const Quad& a, const INT& c);
int (*pos)(const Quad& a);
Quad (*quadgcd)(const Quad& aa, const Quad& bb);
Quad (*quadbezout)(const Quad& aa, const Quad& bb, Quad& xx, Quad& yy);
Quad (*quadconj)(const Quad& a);
Quad mult0(const Quad& a, const Quad& b)
{
return Quad(fmms(a.r,b.r, Quad::n*a.i,b.i), fmma(a.r,b.i, a.i,b.r), a.nm*b.nm);
}
Quad mult1(const Quad& a, const Quad& b)
{
INT x = a.i*b.i;
return Quad(fmms(a.r,b.r, Quad::n,x), fmma(a.r,b.i, a.i,b.r) + x, a.nm*b.nm);
}
// compute a*b+c*d
Quad mma(const Quad& a, const Quad& b, const Quad& c, const Quad& d)
{
INT x = fmma(a.i,b.i, c.i,d.i); // a.i*b.i + c.i*d.i;
INT r = fmma(a.r,b.r, c.r,d.r) - Quad::n*x; // a.r*b.r + c.r*d.r - Quad::n*x;
INT i = fmma(a.r,b.i, a.i,b.r) + fmma(c.r,d.i, c.i,d.r); // a.r*b.i + a.i*b.r + c.r*d.i + c.i*d.r;
if (Quad::t)
i += x;
return Quad(r,i);
}
// compute a*b-c*d
Quad mms(const Quad& a, const Quad& b, const Quad& c, const Quad& d)
{
INT x = fmms(a.i,b.i, c.i,d.i); // a.i*b.i - c.i*d.i;
INT r = fmms(a.r,b.r, c.r,d.r) - Quad::n*x; // a.r*b.r - c.r*d.r - Quad::n*x;
INT i = fmma(a.r,b.i, a.i,b.r) - fmma(c.r,d.i, c.i,d.r); // a.r*b.i + a.i*b.r - (c.r*d.i + c.i*d.r);
if (Quad::t)
i += x;
return Quad(r,i);
}
void Quad::addprod(const Quad& a, const Quad& b) // this +=a*b
{
if (a.nm==0 || b.nm==0) return;
INT c = a.i*b.i;
r += fmms(a.r,b.r, n,c);
i += fmma(a.r,b.i, a.i,b.r);
if (t)
i += c;
setnorm();
}
void Quad::addprod(long a, const Quad& b) // this +=a*b
{
if (a==0 || b.nm==0) return;
r += a*b.r;
i += a*b.i;
setnorm();
}
void Quad::subprod(const Quad& a, const Quad& b) // this -=a*b
{
if (a.nm==0 || b.nm==0) return;
INT c = a.i*b.i;
r -= fmms(a.r,b.r, n,c);
i -= fmma(a.r,b.i, a.i,b.r);
if (t)
i -= c;
setnorm();
}
void Quad::subprod(long a, const Quad& b) // this -=a*b
{
if (a==0 || b.nm==0) return;
r -= a*b.r;
i -= a*b.i;
setnorm();
}
Quad qdivi0(const Quad& a, const INT& c) // c>0, // used when t=0
{
return (c>0?
Quad(rounded_division(a.r,c), rounded_division(a.i,c))
:
Quad(rounded_division(-a.r,-c), rounded_division(-a.i,-c)));
}
Quad qdivi1(const Quad& a, const INT& c) // used when t=1
{
INT ansi = (c>0? rounded_division(a.i,c) : rounded_division(-a.i,-c));
return (c>0?
Quad(rounded_division(2*a.r+a.i-c*ansi,2*c), ansi)
:
Quad(rounded_division(-2*a.r-a.i+c*ansi,-2*c), ansi));
}
// static function (one for the class, not per instance)
void Quad::field(long dd, long max)
{
// Clear these in case this is not the first field run
geometry_initialised = 0;
quadunits.clear();
squareunits.clear();
quadprimes.clear();
prime_disc_factors.clear();
all_disc_factors.clear();
class_group.clear();
class_group_2_torsion.clear();
class_group_2_cotorsion.clear();
class_group_2_torsion_gens.clear();
class_group_2_cotorsion_gens.clear();
Qideal_lists::clear();
SD.clear();
d = squarefree_part(dd);
if (d!=dd)
cout << "Replacing d = " << dd << " with " << d << endl;
is_Euclidean = check_field(d, euclidean_fields);
class_number = 0;
if (check_field(d, class_number_one_fields))
class_number=1;
else
if (check_field(d, class_number_two_fields))
class_number=2;
else
if (check_field(d, class_number_three_fields))
class_number=3;
else
if (check_field(d, class_number_four_fields))
class_number=4;
else
if (check_field(d, class_number_five_fields))
class_number=5;
// else class number is set in fill_class_group()
INT odd(d&1?d:d/2); // odd part of d
switch (d%4)
{
case 1:
{
t=0; absdisc=4*d; disc=-absdisc; n=d;
prime_disc_factors.push_back(INT(-4));
quadconj=&quadconj0;
mult=&mult0; qdivi=&qdivi0;
break;
}
case 2:
{
t=0; absdisc=4*d; disc=-absdisc; n=d;
prime_disc_factors.push_back(INT(odd%4==1 ? -8 : 8));
quadconj=&quadconj0;
mult=&mult0; qdivi=&qdivi0;
break;
}
case 3:
default:
{
t=1; absdisc=d; disc=-d; n=(d+1)/4;
quadconj=&quadconj1;
mult=&mult1; qdivi=&qdivi1;
}
}
// append odd primes divisors of d with sign to be 1 mod 4:
vector<INT> pp = pdivs(odd);
std::for_each(pp.begin(), pp.end(), [] (INT& p) { if (p%4==3) p=-p;});
prime_disc_factors.insert(prime_disc_factors.end(), pp.begin(), pp.end());
INT i0(0), i1(1);
w = Quad(i0, i1, n);
zero = Quad(i0,i0, i0);
one = Quad(i1,i0, i1);
shifts_by_one = {-one-w, -one, -one+w, -w, zero, w, one-w, one, one+w};
switch (d) {
case 1: pos=&pos13; name='i'; nunits=4; fundunit=w; break;
case 2: pos=&pos2; name='t'; nunits=2; fundunit=-one; break;
case 3: pos=&pos13; name='w'; nunits=6; fundunit=w; break;
default: pos=&pos2; name='a'; nunits=2; fundunit=-one;
}
quadgcd=&quadgcd_default;
quadbezout=&quadbezout_default;
quadunits.push_back(one);
quadunits.push_back(fundunit);
for(int i=2; i<nunits; i++)
quadunits.push_back(fundunit*quadunits[i-1]);
for(int i=0; 2*i<nunits; i++)
squareunits.push_back(quadunits[2*i]);
maxnorm=max;
if(class_number==1)
initquadprimes();
Quadprimes::init(max);
fill_class_group();
int n2r = class_group_2_rank; // which was set in fill_class_group()
int nchi = (1<<n2r);
for(int chi_index=0; chi_index<nchi; chi_index++)
{
INT D(1);
for (int i=0; i<n2r; i++)
if (bit(chi_index,i)==1)
D *= prime_disc_factors[i];
if (D>1) // negate positive ones except D=1 itself
D = Quad::disc/D;
all_disc_factors.push_back(D);
}
}
// static function (one for the class, not per instance)
void Quad::displayfield(ostream& s, int info2)
{s<<"Q(sqrt("<<-d<<"))\tdiscriminant = "<<disc;
s<<"\tmin poly("<<name<<") = "<<name<<"^2"; if(t) s<<"-"<<name;
s<<"+"<<n<<".\n";
switch (d) {
case 1: case 2: case 3: case 7: case 11: // Euclidean
cout << "Euclidean" << endl;
break;
case 19: case 43: case 67: case 163: // Non-Euclidean
cout << "Non-Euclidean, class number 1" << endl;
break;
default:
cout << "Class number " << class_number << endl;
cout << "Ideal class group representatives: " << class_group << endl;
if (info2)
{
cout << "2-rank of class group = " << class_group_2_rank << endl;
if (class_group_2_rank>0)
{
cout << " discriminant = ";
for (auto di=prime_disc_factors.begin(); di!=prime_disc_factors.end(); ++di)
{
if (di!=prime_disc_factors.begin()) cout<<"*";
cout<<"("<<(*di)<<")";
}
cout<<endl;
cout << " 2-torsion generators " << class_group_2_torsion_gens
<< ", 2-torsion elements "<< class_group_2_torsion<<endl;
cout << " 2-cotorsion generators " << class_group_2_cotorsion_gens
<< ", 2-cotorsion elements "<< class_group_2_cotorsion<<endl;
}
}
}
if (class_number==1)
s<<nquadprimes<<" primes initialised, max norm = " << maxnorm << endl;
}
int Quad::chi(const INT& p)
{
return (p==2? (d%4==3? (d%8==3? -1: +1): 0): legendre(disc,p));
}
Quad makepos(const Quad& a)
{
Quad ans=a;
while(!pos(ans)) {ans*=fundunit; }
return ans;
}
istream& operator>>(istream& s, Quad& x)
{
s >> x.r >> x.i;
x.setnorm();
return s;
}
void Quad::setnorm()
{
nm = fmma(r,r, n*i,i);
if (t) {nm += r*i;};
// INT N = r*r+t*r*i+n*i*i;
// if (N!=nm) cerr<<"In setnorm(): Quad "<<(*this)<<" with r = "<<r<<", i = "<<i
// <<" and norm "<<N<<" has nm = "<<nm<<endl;
// assert(N==nm);
assert (nm>=0);
}
Quad Quad::operator/ (const Quad& b) const
{
if (!::is_zero(b.i))
return qdivi(mult_conj(*this,b), b.nm);
else
return qdivi(*this,b.r);
}
void Quad::operator/=(const Quad& b)
{
if (!::is_zero(b.i))
*this=qdivi(mult_conj(*this,b), b.nm);
else
*this=qdivi(*this,b.r);
}
// binary index i of [I] in C/C^2 (0 <= i < 2^{2-rank}).
//
// "binary" means that the bits of I give the coordinates w.r.t. the
// basis class_group_2_cotorsion_gens
int Quad::ideal_class_mod_squares(const Qideal& I)
{
return find_ideal_class_mod_squares(I, class_group_2_cotorsion);
}
// image (+1/-1) of I under i'th quadratic character (0 <= i < 2^{2-rank})
int Quad::unramified_character(int i, const Qideal& I)
{
return ((i & ideal_class_mod_squares(I)) ? -1 : +1); // bitwise &
}
// Quad::Quad(const bigcomplex& z)
// {bigfloat x=real(z), y=imag(z);
// if(d>1) y/=sqrt(to_bigfloat(d));
// if(d>2) {x-=y; y*=2.0;}
// Iasb(r,x); Iasb(i,y); //Rounded
// //longify(x, r); longify(y, i); //Rounded
// setnorm();
// }
// Quad::operator bigcomplex() const
// {bigfloat x = to_bigfloat(r), y = to_bigfloat(i);
// if(d>2) {y/=2.0; x+=y;}
// if(d>1) y*=sqrt(to_bigfloat(d));
// return bigcomplex(x,y);
// }
int div(const Quad& a, const Quad& b)
{
if (a.nm==0) return (b.nm==0);
if (b.nm==0) return 1;
if (b.nm%a.nm!=0) return 0;
Quad c = mult_conj(b,a);
return (((c.r)%a.nm)==0) && (((c.i)%a.nm)==0);
}
// as above but return quotient b/a when a|b
int div(const Quad& a, const Quad& b, Quad& quo)
{
if (a.nm==0)
{
quo=0;
return (b.nm==0);
}
if (b.nm==0)
{
quo=0;
return 1;
}
if (b.nm%a.nm!=0)
return 0;
Quad c = mult_conj(b,a);
INT qr, qi, rr, ri;
if ( divrem(c.r, a.nm, qr, rr) && divrem(c.i, a.nm, qi, ri) )
{
quo = Quad(qr,qi);
return 1;
}
return 0;
}
int ndiv(const Quad& a, const Quad& b)
{
return !div(a,b);
}
int val(const Quad& factor, const Quad& number)
{
if ((number.nm==0) || (factor.nm<=1))
{
cout << "Error in val(): factor = "<<factor<< " should not be a unit"<<endl;
exit(1);
}
int e = 0; Quad n=number, f=factor, nf;
while (nf=n/f, f*nf==n) {e++; n=nf;}
return e;
}
vector<Quad> residues(const Quad& a)
{
INT norma = a.norm(), m = gcd(a.re(), a.im());
INT rednorma = (norma/m)/m;
vector<Quad> ans;
for(int j=0; j<m*rednorma; j++)
{
INT J(j);
for(int k=0; k<m; k++)
ans.push_back(Quad(J,INT(k))%a);
}
return ans;
}
vector<Quad> invertible_residues(const Quad& a)
{
vector<Quad> res = residues(a);
vector<Quad> ires(res.size()); // will be shrunk later
auto it = std::copy_if(res.begin(), res.end(), ires.begin(),
[a] (const Quad& r) {return coprime(a,r);});
ires.resize(std::distance(ires.begin(), it));
return ires;
}
Quad::operator string() const
{
ostringstream s;
if (i==0)
s<<r;
else
{
if (r==0)
{
if (i==1) ;
else if (i==-1) s << "-";
else s << i;
}
else
{
s<<r;
if(i>0) s<<"+"; else s<<"-";
if (abs(i)>1) s<<abs(i);
}
s<<Quad::name;
}
return s.str();
}
ostream& operator<<(ostream& s, const Quad& a)
{
s << (string)a;
return s;
}
//Functions for computing quad-primes, initializing the vector<Quad>
//quadprimes. NB all primes are "pos" i.e. normalized w.r.t. units
void factorp0(long p, INT& a, INT& b, const INT& d)
// finds a,b s.t. a^2+d*b^2=0 (mod p)
{ int found=0;
for (int ib=1; !found; ib++)
{
b = ib;
found = isqrt(p - d*b*b, a);
}
}
void factorp1(long p, INT& a, INT& b, const INT& d)
// finds a,b s.t. a^2+a*b+((d+1)/4)*b^2=0 (mod p)
{ int found=0; long fourp = 4*p;
for (int ib=1; !found; ib++)
{
b = ib;
found = isqrt(fourp -d*b*b,a);
}
a=(a-b)/2;
}
vector<Quad> Quad::primes_above(long p, int& sig)
{
INT d(Quad::d), P(p), a, b;
Quad pi, piconj;
vector<Quad> list;
sig = Quad::chi(P);
// cout<<"disc = "<<Quad::disc<<", p="<<p<<", chi(p)="<<sig<<endl;
INT i0(0), i1(1), i2(2), i3(3);
switch (sig) {
case 0: // ramified
pi = (d==1 ? Quad(i1,i1, i2) :
d==2 ? Quad(i0,i1, i2) :
d==3 ? Quad(i1,i1, i3) :
Quad(-i1,i2, d));
list.push_back(pi);
break;
case -1: // inert
pi = Quad(P,i0, P*P);
list.push_back(pi);
break;
case 1: // split
if(Quad::t) factorp1(p,a,b,d); else factorp0(p,a,b,d);
pi = makepos(Quad(a,b, P));
piconj = makepos(quadconj(pi));
// We choose the ordering so the HNFs are [p,c,1], [p,c',1] with c<c'
long c = posmod((a%p)*invmod(b,p),p);
if (2*c<p)
{
list.push_back(pi);
list.push_back(piconj);
}
else
{
list.push_back(piconj);
list.push_back(pi);
}
}
//cout << "primes_above("<<p<<") = " << list << endl;
return list;
}
void Quad::initquadprimes()
{
int sig;
vector<Quad> list, list1, list2;
for (primevar pr; pr.ok()&&pr<maxnorm; pr++)
{ long p=pr;
list = Quad::primes_above(p, sig);
if (sig==-1)
{
if (p*p<=maxnorm)
list2.push_back(list[0]);
}
else
list1.insert(list1.end(), list.begin(), list.end());
}
// Now list1 contains the degree 1 primes and list2 the degree 2
// primes, each ordered by norm; we merge these lists so that they
// are still sorted by norm:
auto alpha=list1.begin(), beta=list2.begin();
while ((alpha!=list1.end()) && (beta!=list2.end()))
if ((alpha->nm)<(beta->nm))
quadprimes.push_back(*alpha++);
else
quadprimes.push_back(*beta++);
while (alpha!=list1.end()) quadprimes.push_back(*alpha++);
while ( beta!=list2.end()) quadprimes.push_back(*beta++);
nquadprimes = quadprimes.size();
}
// Quad::fill_class_group() is implemented in qidloop.cc
Quad primdiv(const Quad& a)
{
INT na = a.norm();
if (na<2) return Quad::zero; // must return something!
for ( const auto& p : quadprimes)
{
if (div(p,a)) return p;
INT np = p.norm();
if (np*np>na) return makepos(a);
}
cout<<"No prime divisor found for "<<a<<" so assuming prime!\n";
return makepos(a);
}
vector<Quad> pdivs(const Quad& aa)
{ Quad a=aa; INT norma = a.norm();
vector<Quad> plist; // will hold prime factors
if (norma<2) return plist;
for ( const auto& p : quadprimes)
{
if (div(p,a))
{
plist.push_back(p);
while (div(p,a, a));
norma = a.norm();
if (norma==1) return plist;
}
else
{
INT normp = p.norm();
if (normp*normp>norma)
{
plist.push_back(makepos(a));
return plist;
}
}
}
//In case of p-factors outside range, assume the cofactor is prime:
if (a.norm()>1)
plist.push_back(makepos(a));
return plist;
}
vector<Quad> posdivs(const Quad& a) // all "positive" divisors (up to units)
{
vector<Quad> plist=pdivs(a);
int nu = 1; int nd=nu;
vector<int> elist;
for ( const auto& p : plist)
{
int e=val(p,a);
elist.push_back(e);
nd*=(1+e);
}
vector<Quad> dlist(nd);
dlist[0]=Quad::one;
nd=nu;
auto pr=plist.begin();
auto ei=elist.begin();
for ( ; pr!=plist.end(); ++pr, ++ei)
{
Quad p = *pr;
int e = *ei;
for (int j=0; j<e; j++)
for (int k=0; k<nd; k++)
dlist[nd*(j+1)+k] = makepos(p*dlist[nd*j+k]);
nd*=(e+1);
}
return dlist;
}
vector<Quad> alldivs(const Quad& a) // all divisors
{
vector<Quad> plist=pdivs(a);
int nu = Quad::nunits; int nd=nu;
vector<int> elist;
for ( const auto& p :plist)
{
int e=val(p,a);
elist.push_back(e);
nd*=(1+e);
}
vector<Quad> dlist(nd);
for(int i=0; i<nu; i++) dlist[i]=quadunits[i];
nd=nu;
auto pr=plist.begin();
auto ei=elist.begin();
for ( ; pr!=plist.end(); ++pr, ++ei)
{
Quad p = *pr;
int e = *ei;
for (int j=0; j<e; j++)
for (int k=0; k<nd; k++)
dlist[nd*(j+1)+k] = p*dlist[nd*j+k];
nd*=(e+1);
}
return dlist;
}
vector<Quad> sqdivs(const Quad& a) // all divisors whose square divides a, up to +/-
{
vector<Quad> plist=pdivs(a);
int nu = Quad::nunits/2; int nd=nu;
vector<int> elist;
for ( const auto& p : plist)
{
int e=val(p,a)/2;
elist.push_back(e);
nd*=(1+e);
}
vector<Quad> dlist(nd);
for(int i=0; i<nu; i++) dlist[i]=quadunits[i];
nd=nu;
auto pr=plist.begin();
auto ei=elist.begin();
for ( ; pr!=plist.end(); ++pr, ++ei)
{
Quad p = *pr;
int e = *ei;
for (int j=0; j<e; j++)
for (int k=0; k<nd; k++)
dlist[nd*(j+1)+k] = p*dlist[nd*j+k];
nd*=(e+1);
}
return dlist;
}
vector<Quad> sqfreedivs(const Quad& a) // all square-free divisors
{
vector<Quad> plist=pdivs(a);
int nu = 2; int nd=pow(2,plist.size()+1);
vector<int> elist(plist.size(), 1);
vector<Quad> dlist(nd);
for(int i=0; i<nu; i++) dlist[i]=quadunits[i];
nd=nu;
auto pr=plist.begin();
auto ei=elist.begin();
for ( ; pr!=plist.end(); ++pr, ++ei)
{
Quad p = *pr;
int e = *ei;
for (int j=0; j<e; j++)
for (int k=0; k<nd; k++)
dlist[nd*(j+1)+k] = p*dlist[nd*j+k];
nd*=(e+1);
}
return dlist;
}
Quad invmod(const Quad& a, const Quad& p)
{Quad x,y;
Quad g=quadbezout(a,p,x,y);
if (g==Quad::one) return x;
else {
cerr<<"invmod called with "<<a<<" and "<<p<<" -- not coprime!"<<endl;
return Quad::zero;
}
}
int invertible(const Quad& a, const Quad& b, Quad& inverse)
{ Quad y; Quad g = quadbezout(a,b,inverse,y);
return g.is_one();
}
//#define test_reduce
// returns n = round(true_real_part_of(alpha/beta)), so alpha-n*beta
// is reduced mod Z<beta>
INT nearest_INT_to_Quad_quotient ( const Quad& alpha, const Quad& beta)
{
return rounded_division(racb(alpha,beta), beta.norm());
}
// reduction of gamma modulo Z<alpha,beta>
Quad reduce_mod_zbasis(const Quad& gamma, const Quad& alpha, const Quad& beta)
{
#ifdef test_reduce
cout << "reduction of "<<gamma<<" mod <"<<alpha<<","<<beta<<">"<< flush;
#endif
INT d = iacb(beta, alpha); // = (beta*quadconj(alpha)).im();
assert (d>0);
INT x = rounded_division(-iacb(gamma, beta), d);
INT y = rounded_division( iacb(gamma, alpha), d);
Quad ans = gamma - (x*alpha + y*beta);
#ifdef test_reduce
cout << " is "<< ans << " (d="<<d<<", x="<<x<<", y="<<y<<")"<<endl;
cout << " x*alpha = "<< x*alpha << endl;
cout << " y*beta = "<< y*beta << endl;
cout << " x*alpha+y*beta = "<< x*alpha+y*beta << endl;
#endif
INT gn = ans.norm();
vector<Quad> tests = {ans+alpha, ans-alpha, ans+beta, ans-beta,
ans-alpha-beta, ans-alpha+beta, ans+alpha-beta, ans+alpha+beta};
for ( const auto& t : tests)
{
if (gn > t.norm())
{
#ifdef test_reduce
cout<<"reduction "<<ans<<" has larger norm than shift "<<t<<", switching"<<endl;
#endif
ans = t;
gn = t.norm();
}
}
return ans;
}
#undef test_reduce
// Replace alpha, beta by an SL2Z-equivalent pair with the same Z-span.
// The new alpha has the smallest norm.
void sl2z_reduce(Quad& alpha, Quad& beta)
{
#ifdef test_reduce
cout<<"SL2Z-reducing ["<<alpha<<","<<beta<<"]..."<<endl;
#endif
int s=1;
while (s)
{
s = 0; // will be set to 1 if anything changes
INT n = nearest_INT_to_Quad_quotient(beta,alpha);
if(n!=0)
{
s=1;
beta -= n*alpha;
#ifdef test_reduce
cout<<" -- shift by " << n <<": alpha="<<alpha<<", beta="<<beta<< endl;
#endif
}
if (beta.nm < alpha.nm)
{
s=1;
Quad temp = -alpha; alpha = beta; beta = temp;
#ifdef test_reduce
cout<<" -- invert: alpha="<<alpha<<", beta="<<beta<< endl;
#endif
}
}
// alpha is now the non-zero Quad of least norm in the lattice [alpha,beta]
// and beta/alpha is in the closed fundamental region
// We want to orient so that beta/alpha has positive imaginary part
// e.g. so that the basis [1,w] is reduced
// NB im(x+y*w)>0 iff y>0 for both t=0 and t=1 cases
if (iacb(beta, alpha) < 0)
beta=-beta;
// If we wanted to make sure that beta/alpha is unambiguously defined when on the boundary:
// if (2*racm(beta, alpha) == -alpha.norm())
// beta+=alpha;
#ifdef test_reduce
cout<<"After reduction, alpha="<<alpha<<", beta="<<beta
<<" with norms "<<alpha.norm()<<", "<<beta.norm()<<endl;
assert (alpha.norm()<=beta.norm());
//assert (nearest_INT_to_Quad_quotient(alpha,beta)==0);
assert (iacb(beta, alpha) > 0);
#endif
}
// Replace alpha, beta by an SL2Z-equivalent pair with the same Z-span.
// The new alpha has the smallest norm.
// U holds the unimodular transform (not implemented for FLINT INTs as the unimod class uses NTL ZZ)
// void sl2z_reduce(Quad& alpha, Quad& beta, unimod&U)
// {
// Quad alpha0=alpha, beta0=beta;
// INT U11, U12, U21, U22;
// #ifdef test_reduce
// cout<<"SL2Z-reducing ["<<alpha<<","<<beta<<"]..."<<endl;
// #endif
// U.reset(); // to the identity
// int s=1;
// while (s)
// {
// s = 0; // will be set to 1 if anything changes
// INT n = nearest_INT_to_Quad_quotient(beta,alpha);
// if(n!=0)
// {
// s=1;
// beta -= n*alpha;
// U.y_shift(INT(n));
// #ifdef test_reduce
// cout<<" -- shift by " << n <<": alpha="<<alpha<<", beta="<<beta<< endl;
// #endif
// U11 = I2long(U(1,1)); U12 = I2long(U(1,2));
// U21 = I2long(U(2,1)); U22 = I2long(U(2,2));
// assert (U11*alpha+U12*beta == alpha0);
// assert (U21*alpha+U22*beta == beta0);
// }
// if (beta.nm < alpha.nm)
// {
// s=1;
// Quad temp = -alpha; alpha = beta; beta = temp;
// U.invert();
// #ifdef test_reduce
// cout<<" -- invert: alpha="<<alpha<<", beta="<<beta<< endl;
// U11 = I2long(U(1,1)); U12 = I2long(U(1,2));
// U21 = I2long(U(2,1)); U22 = I2long(U(2,2));
// assert (U11*alpha+U12*beta == alpha0);
// assert (U21*alpha+U22*beta == beta0);
// #endif
// }
// }
// // alpha is now the non-zero Quad of least norm in the lattice [alpha,beta]
// // and beta/alpha is in the closed fundamental region
// // We want to orient so that beta/alpha has positive imaginary part
// // e.g. so that the basis [1,w] is reduced
// // NB im(x+y*w)>0 iff y>0 for both t=0 and t=1 cases
// if (iacb(beta,alpha) < 0)
// beta=-beta;
// // If we wanted to make sure that beta/alpha is unambiguously defined when on the boundary:
// // if (2*tacm(beta, alpha) == -alpha.norm())
// // beta+=alpha;
// #ifdef test_reduce
// cout<<"After reduction by U="<<U<<", alpha="<<alpha<<", beta="<<beta
// <<" with norms "<<alpha.norm()<<", "<<beta.norm()<<endl;
// U11 = I2long(U(1,1)); U12 = I2long(U(1,2));
// U21 = I2long(U(2,1)); U22 = I2long(U(2,2));
// assert (U11*alpha+U12*beta == alpha0);
// assert (U21*alpha+U22*beta == beta0);
// assert (alpha.norm()<=beta.norm());
// assert (nearest_INT_to_Quad_quotient(alpha,beta)==0);
// assert (iacb(beta, alpha) > 0);
// #endif
// }
// HNF of ideal (alpha) as a triple [a c d] where [a,c+d*w] is a Z-basis with
//
// a,d>0; c>=0
// N=a*d = Norm(alpha)