Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/matrix/benchmark.py
4049 views
1
"""
2
Benchmarks for matrices
3
4
This file has many functions for computing timing benchmarks
5
of various methods for random matrices with given bounds for
6
the entries. The systems supported are Sage and Magma.
7
8
The basic command syntax is as follows::
9
10
sage: import sage.matrix.benchmark as b
11
sage: b.report([b.det_ZZ], 'Test', systems=['sage'])
12
======================================================================
13
Test
14
======================================================================
15
...
16
======================================================================
17
"""
18
from sage.all import *
19
20
verbose = False
21
22
timeout = 60
23
24
def report(F, title, systems = ['sage', 'magma'], **kwds):
25
"""
26
Run benchmarks with default arguments for each function in the list F.
27
28
INPUT:
29
30
- ``F`` - a list of callables used for benchmarking
31
- ``title`` - a string describing this report
32
- ``systems`` - a list of systems (supported entries are 'sage' and 'magma')
33
- ``**kwds`` - keyword arguments passed to all functions in ``F``
34
35
EXAMPLES::
36
37
sage: import sage.matrix.benchmark as b
38
sage: b.report([b.det_ZZ], 'Test', systems=['sage'])
39
======================================================================
40
Test
41
======================================================================
42
...
43
======================================================================
44
"""
45
if len(systems) > 2:
46
raise NotImplementedError, "at most two systems ('sage' or 'magma')"
47
print '='*70
48
print ' '*10 + title
49
print '='*70
50
os.system('uname -a')
51
print '\n'
52
for f in F:
53
print "-"*70
54
print f.__doc__.strip()
55
print ('%15s'*len(systems))%tuple(systems)
56
w = []
57
for s in systems:
58
alarm(timeout)
59
try:
60
t = f(system=s, **kwds)
61
except KeyboardInterrupt:
62
t = -timeout
63
alarm(0)
64
w.append(float(t))
65
if len(w) > 1:
66
if w[1] == 0:
67
w.append(0.0)
68
else:
69
w.append(w[0]/w[1])
70
71
w = tuple(w)
72
print ('%15.3f'*len(w))%w
73
print '='*70
74
75
76
#######################################################################
77
# Dense Benchmarks over ZZ
78
#######################################################################
79
80
def report_ZZ(**kwds):
81
"""
82
Reports all the benchmarks for integer matrices and few
83
rational matrices.
84
85
INPUT:
86
87
- ``**kwds`` - passed through to :func:`report`
88
89
EXAMPLES::
90
91
sage: import sage.matrix.benchmark as b
92
sage: b.report_ZZ(systems=['sage']) # long time (15s on sage.math, 2012)
93
======================================================================
94
Dense benchmarks over ZZ
95
======================================================================
96
...
97
======================================================================
98
"""
99
F = [vecmat_ZZ, rank_ZZ, rank2_ZZ, charpoly_ZZ, smithform_ZZ,
100
det_ZZ, det_QQ, matrix_multiply_ZZ, matrix_add_ZZ,
101
matrix_add_ZZ_2,
102
nullspace_ZZ]
103
104
title = 'Dense benchmarks over ZZ'
105
report(F, title, **kwds)
106
107
# Integer Nullspace
108
109
def nullspace_ZZ(n=200, min=0, max=2**32, system='sage'):
110
"""
111
Nullspace over ZZ:
112
Given a n+1 x n matrix over ZZ with random entries
113
between min and max, compute the nullspace.
114
115
INPUT:
116
117
- ``n`` - matrix dimension (default: ``200``)
118
- ``min`` - minimal value for entries of matrix (default: ``0``)
119
- ``max`` - maximal value for entries of matrix (default: ``2**32``)
120
- ``system`` - either 'sage' or 'magma' (default: 'sage')
121
122
EXAMPLES::
123
124
sage: import sage.matrix.benchmark as b
125
sage: ts = b.nullspace_ZZ(200)
126
sage: tm = b.nullspace_ZZ(200, system='magma') # optional - magma
127
"""
128
if system == 'sage':
129
A = random_matrix(ZZ, n+1, n, x=min, y=max+1).change_ring(QQ)
130
t = cputime()
131
v = A.kernel()
132
return cputime(t)
133
elif system == 'magma':
134
code = """
135
n := %s;
136
A := RMatrixSpace(RationalField(), n+1,n)![Random(%s,%s) : i in [1..n*(n+1)]];
137
t := Cputime();
138
K := Kernel(A);
139
s := Cputime(t);
140
"""%(n,min,max)
141
if verbose: print code
142
magma.eval(code)
143
return float(magma.eval('s'))
144
else:
145
raise ValueError, 'unknown system "%s"'%system
146
147
148
def charpoly_ZZ(n=100, min=0, max=9, system='sage'):
149
"""
150
Characteristic polynomial over ZZ:
151
Given a n x n matrix over ZZ with random entries between min and
152
max, compute the charpoly.
153
154
INPUT:
155
156
- ``n`` - matrix dimension (default: ``100``)
157
- ``min`` - minimal value for entries of matrix (default: ``0``)
158
- ``max`` - maximal value for entries of matrix (default: ``9``)
159
- ``system`` - either 'sage' or 'magma' (default: 'sage')
160
161
EXAMPLES::
162
163
sage: import sage.matrix.benchmark as b
164
sage: ts = b.charpoly_ZZ(100)
165
sage: tm = b.charpoly_ZZ(100, system='magma') # optional - magma
166
"""
167
if system == 'sage':
168
A = random_matrix(ZZ, n, n, x=min, y=max+1)
169
t = cputime()
170
v = A.charpoly()
171
return cputime(t)
172
elif system == 'magma':
173
code = """
174
n := %s;
175
A := MatrixAlgebra(IntegerRing(), n)![Random(%s,%s) : i in [1..n^2]];
176
t := Cputime();
177
K := CharacteristicPolynomial(A);
178
s := Cputime(t);
179
"""%(n,min,max)
180
if verbose: print code
181
magma.eval(code)
182
return float(magma.eval('s'))
183
else:
184
raise ValueError, 'unknown system "%s"'%system
185
186
187
def rank_ZZ(n=700, min=0, max=9, system='sage'):
188
"""
189
Rank over ZZ:
190
Given a n x (n+10) matrix over ZZ with random entries
191
between min and max, compute the rank.
192
193
INPUT:
194
195
- ``n`` - matrix dimension (default: ``700``)
196
- ``min`` - minimal value for entries of matrix (default: ``0``)
197
- ``max`` - maximal value for entries of matrix (default: ``9``)
198
- ``system`` - either 'sage' or 'magma' (default: 'sage')
199
200
EXAMPLES::
201
202
sage: import sage.matrix.benchmark as b
203
sage: ts = b.rank_ZZ(300)
204
sage: tm = b.rank_ZZ(300, system='magma') # optional - magma
205
"""
206
if system == 'sage':
207
A = random_matrix(ZZ, n, n+10, x=min, y=max+1)
208
t = cputime()
209
v = A.rank()
210
return cputime(t)
211
elif system == 'magma':
212
code = """
213
n := %s;
214
A := RMatrixSpace(IntegerRing(), n, n+10)![Random(%s,%s) : i in [1..n*(n+10)]];
215
t := Cputime();
216
K := Rank(A);
217
s := Cputime(t);
218
"""%(n,min,max)
219
if verbose: print code
220
magma.eval(code)
221
return float(magma.eval('s'))
222
else:
223
raise ValueError, 'unknown system "%s"'%system
224
225
def rank2_ZZ(n=400, min=0, max=2**64, system='sage'):
226
"""
227
Rank 2 over ZZ:
228
Given a (n + 10) x n matrix over ZZ with random entries
229
between min and max, compute the rank.
230
231
INPUT:
232
233
- ``n`` - matrix dimension (default: ``400``)
234
- ``min`` - minimal value for entries of matrix (default: ``0``)
235
- ``max`` - maximal value for entries of matrix (default: ``2**64``)
236
- ``system`` - either 'sage' or 'magma' (default: 'sage')
237
238
EXAMPLES::
239
240
sage: import sage.matrix.benchmark as b
241
sage: ts = b.rank2_ZZ(300)
242
sage: tm = b.rank2_ZZ(300, system='magma') # optional - magma
243
"""
244
if system == 'sage':
245
A = random_matrix(ZZ, n+10, n, x=min, y=max+1)
246
t = cputime()
247
v = A.rank()
248
return cputime(t)
249
elif system == 'magma':
250
code = """
251
n := %s;
252
A := RMatrixSpace(IntegerRing(), n+10, n)![Random(%s,%s) : i in [1..n*(n+10)]];
253
t := Cputime();
254
K := Rank(A);
255
s := Cputime(t);
256
"""%(n,min,max)
257
if verbose: print code
258
magma.eval(code)
259
return float(magma.eval('s'))
260
else:
261
raise ValueError, 'unknown system "%s"'%system
262
263
# Smith Form
264
265
def smithform_ZZ(n=128, min=0, max=9, system='sage'):
266
"""
267
Smith Form over ZZ:
268
Given a n x n matrix over ZZ with random entries
269
between min and max, compute the Smith normal form.
270
271
INPUT:
272
273
- ``n`` - matrix dimension (default: ``128``)
274
- ``min`` - minimal value for entries of matrix (default: ``0``)
275
- ``max`` - maximal value for entries of matrix (default: ``9``)
276
- ``system`` - either 'sage' or 'magma' (default: 'sage')
277
278
EXAMPLES::
279
280
sage: import sage.matrix.benchmark as b
281
sage: ts = b.smithform_ZZ(100)
282
sage: tm = b.smithform_ZZ(100, system='magma') # optional - magma
283
"""
284
if system == 'sage':
285
A = random_matrix(ZZ, n, n, x=min, y=max+1)
286
t = cputime()
287
v = A.elementary_divisors()
288
return cputime(t)
289
elif system == 'magma':
290
code = """
291
n := %s;
292
A := MatrixAlgebra(IntegerRing(), n)![Random(%s,%s) : i in [1..n^2]];
293
t := Cputime();
294
K := ElementaryDivisors(A);
295
s := Cputime(t);
296
"""%(n,min,max)
297
if verbose: print code
298
magma.eval(code)
299
return float(magma.eval('s'))
300
else:
301
raise ValueError, 'unknown system "%s"'%system
302
303
304
def matrix_multiply_ZZ(n=300, min=-9, max=9, system='sage', times=1):
305
"""
306
Matrix multiplication over ZZ
307
Given an n x n matrix A over ZZ with random entries
308
between min and max, inclusive, compute A * (A+1).
309
310
INPUT:
311
312
- ``n`` - matrix dimension (default: ``300``)
313
- ``min`` - minimal value for entries of matrix (default: ``-9``)
314
- ``max`` - maximal value for entries of matrix (default: ``9``)
315
- ``system`` - either 'sage' or 'magma' (default: 'sage')
316
- ``times`` - number of experiments (default: ``1``)
317
318
EXAMPLES::
319
320
sage: import sage.matrix.benchmark as b
321
sage: ts = b.matrix_multiply_ZZ(200)
322
sage: tm = b.matrix_multiply_ZZ(200, system='magma') # optional - magma
323
"""
324
if system == 'sage':
325
A = random_matrix(ZZ, n, n, x=min, y=max+1)
326
B = A + 1
327
t = cputime()
328
for z in range(times):
329
v = A * B
330
return cputime(t)/times
331
elif system == 'magma':
332
code = """
333
n := %s;
334
A := MatrixAlgebra(IntegerRing(), n)![Random(%s,%s) : i in [1..n^2]];
335
B := A + 1;
336
t := Cputime();
337
for z in [1..%s] do
338
K := A * B;
339
end for;
340
s := Cputime(t);
341
"""%(n,min,max,times)
342
if verbose: print code
343
magma.eval(code)
344
return float(magma.eval('s'))/times
345
else:
346
raise ValueError, 'unknown system "%s"'%system
347
348
def matrix_add_ZZ(n=200, min=-9, max=9, system='sage', times=50):
349
"""
350
Matrix addition over ZZ
351
Given an n x n matrix A and B over ZZ with random entries between
352
``min`` and ``max``, inclusive, compute A + B ``times`` times.
353
354
INPUT:
355
356
- ``n`` - matrix dimension (default: ``200``)
357
- ``min`` - minimal value for entries of matrix (default: ``-9``)
358
- ``max`` - maximal value for entries of matrix (default: ``9``)
359
- ``system`` - either 'sage' or 'magma' (default: 'sage')
360
- ``times`` - number of experiments (default: ``50``)
361
362
EXAMPLES::
363
364
sage: import sage.matrix.benchmark as b
365
sage: ts = b.matrix_add_ZZ(200)
366
sage: tm = b.matrix_add_ZZ(200, system='magma') # optional - magma
367
"""
368
if system == 'sage':
369
A = random_matrix(ZZ, n, n, x=min, y=max+1)
370
B = random_matrix(ZZ, n, n, x=min, y=max+1)
371
t = cputime()
372
for z in range(times):
373
v = A + B
374
return cputime(t)/times
375
elif system == 'magma':
376
code = """
377
n := %s;
378
min := %s;
379
max := %s;
380
A := MatrixAlgebra(IntegerRing(), n)![Random(min,max) : i in [1..n^2]];
381
B := MatrixAlgebra(IntegerRing(), n)![Random(min,max) : i in [1..n^2]];
382
t := Cputime();
383
for z in [1..%s] do
384
K := A + B;
385
end for;
386
s := Cputime(t);
387
"""%(n,min,max,times)
388
if verbose: print code
389
magma.eval(code)
390
return float(magma.eval('s'))/times
391
else:
392
raise ValueError, 'unknown system "%s"'%system
393
394
def matrix_add_ZZ_2(n=200, bits=16, system='sage', times=50):
395
"""
396
Matrix addition over ZZ.
397
Given an n x n matrix A and B over ZZ with random ``bits``-bit
398
entries, compute A + B.
399
400
INPUT:
401
402
- ``n`` - matrix dimension (default: ``200``)
403
- ``bits`` - bitsize of entries
404
- ``system`` - either 'sage' or 'magma' (default: 'sage')
405
- ``times`` - number of experiments (default: ``50``)
406
407
EXAMPLES::
408
409
sage: import sage.matrix.benchmark as b
410
sage: ts = b.matrix_add_ZZ_2(200)
411
sage: tm = b.matrix_add_ZZ_2(200, system='magma') # optional - magma
412
"""
413
b = 2**bits
414
return matrix_add_ZZ(n=n, min=-b, max=b,system=system, times=times)
415
416
def det_ZZ(n=200, min=1, max=100, system='sage'):
417
"""
418
Dense integer determinant over ZZ.
419
Given an n x n matrix A over ZZ with random entries
420
between min and max, inclusive, compute det(A).
421
422
INPUT:
423
424
- ``n`` - matrix dimension (default: ``200``)
425
- ``min`` - minimal value for entries of matrix (default: ``1``)
426
- ``max`` - maximal value for entries of matrix (default: ``100``)
427
- ``system`` - either 'sage' or 'magma' (default: 'sage')
428
429
EXAMPLES::
430
431
sage: import sage.matrix.benchmark as b
432
sage: ts = b.det_ZZ(200)
433
sage: tm = b.det_ZZ(200, system='magma') # optional - magma
434
"""
435
if system == 'sage':
436
A = random_matrix(ZZ, n, n, x=min, y=max+1)
437
t = cputime()
438
d = A.determinant()
439
return cputime(t)
440
elif system == 'magma':
441
code = """
442
n := %s;
443
A := MatrixAlgebra(IntegerRing(), n)![Random(%s,%s) : i in [1..n^2]];
444
t := Cputime();
445
d := Determinant(A);
446
s := Cputime(t);
447
"""%(n,min,max)
448
if verbose: print code
449
magma.eval(code)
450
return float(magma.eval('s'))
451
else:
452
raise ValueError, 'unknown system "%s"'%system
453
454
455
def det_QQ(n=300, num_bound=10, den_bound=10, system='sage'):
456
"""
457
Dense rational determinant over QQ.
458
Given an n x n matrix A over QQ with random entries
459
with numerator bound and denominator bound, compute det(A).
460
461
INPUT:
462
463
- ``n`` - matrix dimension (default: ``200``)
464
- ``num_bound`` - numerator bound, inclusive (default: ``10``)
465
- ``den_bound`` - denominator bound, inclusive (default: ``10``)
466
- ``system`` - either 'sage' or 'magma' (default: 'sage')
467
468
EXAMPLES::
469
470
sage: import sage.matrix.benchmark as b
471
sage: ts = b.det_QQ(200)
472
sage: ts = b.det_QQ(10, num_bound=100000, den_bound=10000)
473
sage: tm = b.det_QQ(200, system='magma') # optional - magma
474
"""
475
if system == 'sage':
476
A = random_matrix(QQ, n, n, num_bound=num_bound, den_bound=den_bound)
477
t = cputime()
478
d = A.determinant()
479
return cputime(t)
480
elif system == 'magma':
481
code = """
482
n := %s;
483
A := MatrixAlgebra(RationalField(), n)![Random(%s,%s)/Random(1,%s) : i in [1..n^2]];
484
t := Cputime();
485
d := Determinant(A);
486
s := Cputime(t);
487
"""%(n,-num_bound, num_bound, den_bound)
488
if verbose: print code
489
magma.eval(code)
490
return float(magma.eval('s'))
491
else:
492
raise ValueError, 'unknown system "%s"'%system
493
494
495
def vecmat_ZZ(n=300, min=-9, max=9, system='sage', times=200):
496
"""
497
Vector matrix multiplication over ZZ.
498
499
Given an n x n matrix A over ZZ with random entries
500
between min and max, inclusive, and v the first row of A,
501
compute the product v * A.
502
503
INPUT:
504
505
- ``n`` - matrix dimension (default: ``300``)
506
- ``min`` - minimal value for entries of matrix (default: ``-9``)
507
- ``max`` - maximal value for entries of matrix (default: ``9``)
508
- ``system`` - either 'sage' or 'magma' (default: 'sage')
509
- ``times`` - number of runs (default: ``200``)
510
511
EXAMPLES::
512
513
sage: import sage.matrix.benchmark as b
514
sage: ts = b.vecmat_ZZ(300)
515
sage: tm = b.vecmat_ZZ(300, system='magma') # optional - magma
516
"""
517
if system == 'sage':
518
A = random_matrix(ZZ, n, n, x=min, y=max+1)
519
v = A.row(0)
520
t = cputime()
521
for z in range(times):
522
w = v * A
523
return cputime(t)/times
524
elif system == 'magma':
525
code = """
526
n := %s;
527
A := MatrixAlgebra(IntegerRing(), n)![Random(%s,%s) : i in [1..n^2]];
528
v := A[1];
529
t := Cputime();
530
for z in [1..%s] do
531
K := v * A;
532
end for;
533
s := Cputime(t);
534
"""%(n,min,max,times)
535
if verbose: print code
536
magma.eval(code)
537
return float(magma.eval('s'))/times
538
else:
539
raise ValueError, 'unknown system "%s"'%system
540
541
542
543
#######################################################################
544
# Dense Benchmarks over GF(p), for small p.
545
#######################################################################
546
547
def report_GF(p=16411, **kwds):
548
"""
549
Runs all the reports for finite field matrix operations, for
550
prime p=16411.
551
552
INPUT:
553
554
- ``p`` - ignored
555
- ``**kwds`` - passed through to :func:`report`
556
557
.. note::
558
559
right now, even though p is an input, it is being ignored! If
560
you need to check the performance for other primes, you can
561
call individual benchmark functions.
562
563
EXAMPLES::
564
565
sage: import sage.matrix.benchmark as b
566
sage: b.report_GF(systems=['sage'])
567
======================================================================
568
Dense benchmarks over GF with prime 16411
569
======================================================================
570
...
571
======================================================================
572
"""
573
F = [rank_GF, rank2_GF, nullspace_GF, charpoly_GF,
574
matrix_multiply_GF, det_GF]
575
title = 'Dense benchmarks over GF with prime %i' % p
576
report(F, title, **kwds)
577
578
# Nullspace over GF
579
580
def nullspace_GF(n=300, p=16411, system='sage'):
581
"""
582
Given a n+1 x n matrix over GF(p) with random
583
entries, compute the nullspace.
584
585
INPUT:
586
587
- ``n`` - matrix dimension (default: 300)
588
- ``p`` - prime number (default: ``16411``)
589
- ``system`` - either 'magma' or 'sage' (default: 'sage')
590
591
EXAMPLES::
592
593
sage: import sage.matrix.benchmark as b
594
sage: ts = b.nullspace_GF(300)
595
sage: tm = b.nullspace_GF(300, system='magma') # optional - magma
596
"""
597
if system == 'sage':
598
A = random_matrix(GF(p), n, n+1)
599
t = cputime()
600
v = A.kernel()
601
return cputime(t)
602
elif system == 'magma':
603
code = """
604
n := %s;
605
A := Random(RMatrixSpace(GF(%s), n, n+1));
606
t := Cputime();
607
K := Kernel(A);
608
s := Cputime(t);
609
"""%(n,p)
610
if verbose: print code
611
magma.eval(code)
612
return magma.eval('s')
613
else:
614
raise ValueError, 'unknown system "%s"'%system
615
616
617
# Characteristic Polynomial over GF
618
619
def charpoly_GF(n=100, p=16411, system='sage'):
620
"""
621
Given a n x n matrix over GF with random entries, compute the
622
charpoly.
623
624
INPUT:
625
626
- ``n`` - matrix dimension (default: 100)
627
- ``p`` - prime number (default: ``16411``)
628
- ``system`` - either 'magma' or 'sage' (default: 'sage')
629
630
EXAMPLES::
631
632
sage: import sage.matrix.benchmark as b
633
sage: ts = b.charpoly_GF(100)
634
sage: tm = b.charpoly_GF(100, system='magma') # optional - magma
635
"""
636
if system == 'sage':
637
A = random_matrix(GF(p), n, n)
638
t = cputime()
639
v = A.charpoly()
640
return cputime(t)
641
elif system == 'magma':
642
code = """
643
n := %s;
644
A := Random(MatrixAlgebra(GF(%s), n));
645
t := Cputime();
646
K := CharacteristicPolynomial(A);
647
s := Cputime(t);
648
"""%(n,p)
649
if verbose: print code
650
magma.eval(code)
651
return magma.eval('s')
652
else:
653
raise ValueError, 'unknown system "%s"'%system
654
655
def matrix_add_GF(n=1000, p=16411, system='sage',times=100):
656
"""
657
Given two n x n matrix over GF(p) with random entries, add them.
658
659
INPUT:
660
661
- ``n`` - matrix dimension (default: 300)
662
- ``p`` - prime number (default: ``16411``)
663
- ``system`` - either 'magma' or 'sage' (default: 'sage')
664
- ``times`` - number of experiments (default: ``100``)
665
666
EXAMPLES::
667
668
sage: import sage.matrix.benchmark as b
669
sage: ts = b.matrix_add_GF(500, p=19)
670
sage: tm = b.matrix_add_GF(500, p=19, system='magma') # optional - magma
671
"""
672
if system == 'sage':
673
A = random_matrix(GF(p), n, n)
674
B = random_matrix(GF(p), n, n)
675
t = cputime()
676
for n in range(times):
677
v = A + B
678
return cputime(t)
679
elif system == 'magma':
680
code = """
681
n := %s;
682
A := Random(MatrixAlgebra(GF(%s), n));
683
B := Random(MatrixAlgebra(GF(%s), n));
684
t := Cputime();
685
for z in [1..%s] do
686
K := A + B;
687
end for;
688
s := Cputime(t);
689
"""%(n,p,p,times)
690
if verbose: print code
691
magma.eval(code)
692
return magma.eval('s')
693
else:
694
raise ValueError, 'unknown system "%s"'%system
695
696
697
698
# Matrix multiplication over GF(p)
699
700
def matrix_multiply_GF(n=100, p=16411, system='sage', times=3):
701
"""
702
Given an n x n matrix A over GF(p) with random entries, compute
703
A * (A+1).
704
705
INPUT:
706
707
- ``n`` - matrix dimension (default: 100)
708
- ``p`` - prime number (default: ``16411``)
709
- ``system`` - either 'magma' or 'sage' (default: 'sage')
710
- ``times`` - number of experiments (default: ``3``)
711
712
EXAMPLES::
713
714
sage: import sage.matrix.benchmark as b
715
sage: ts = b.matrix_multiply_GF(100, p=19)
716
sage: tm = b.matrix_multiply_GF(100, p=19, system='magma') # optional - magma
717
"""
718
if system == 'sage':
719
A = random_matrix(GF(p), n)
720
B = A + 1
721
t = cputime()
722
for n in range(times):
723
v = A * B
724
return cputime(t) / times
725
elif system == 'magma':
726
code = """
727
n := %s;
728
A := Random(MatrixAlgebra(GF(%s), n));
729
B := A + 1;
730
t := Cputime();
731
for z in [1..%s] do
732
K := A * B;
733
end for;
734
s := Cputime(t);
735
"""%(n,p,times)
736
if verbose: print code
737
magma.eval(code)
738
return float(magma.eval('s'))/times
739
else:
740
raise ValueError, 'unknown system "%s"'%system
741
742
743
def rank_GF(n=500, p=16411, system='sage'):
744
"""
745
Rank over GF(p):
746
Given a n x (n+10) matrix over GF(p) with random entries, compute the rank.
747
748
INPUT:
749
750
- ``n`` - matrix dimension (default: 300)
751
- ``p`` - prime number (default: ``16411``)
752
- ``system`` - either 'magma' or 'sage' (default: 'sage')
753
754
EXAMPLES::
755
756
sage: import sage.matrix.benchmark as b
757
sage: ts = b.rank_GF(1000)
758
sage: tm = b.rank_GF(1000, system='magma') # optional - magma
759
"""
760
if system == 'sage':
761
A = random_matrix(GF(p), n, n+10)
762
t = cputime()
763
v = A.rank()
764
return cputime(t)
765
elif system == 'magma':
766
code = """
767
n := %s;
768
A := Random(MatrixAlgebra(GF(%s), n));
769
t := Cputime();
770
K := Rank(A);
771
s := Cputime(t);
772
"""%(n,p)
773
if verbose: print code
774
magma.eval(code)
775
return float(magma.eval('s'))
776
else:
777
raise ValueError, 'unknown system "%s"'%system
778
779
def rank2_GF(n=500, p=16411, system='sage'):
780
"""
781
Rank over GF(p): Given a (n + 10) x n matrix over GF(p) with
782
random entries, compute the rank.
783
784
INPUT:
785
786
- ``n`` - matrix dimension (default: 300)
787
- ``p`` - prime number (default: ``16411``)
788
- ``system`` - either 'magma' or 'sage' (default: 'sage')
789
790
EXAMPLES::
791
792
sage: import sage.matrix.benchmark as b
793
sage: ts = b.rank2_GF(500)
794
sage: tm = b.rank2_GF(500, system='magma') # optional - magma
795
"""
796
if system == 'sage':
797
A = random_matrix(GF(p), n+10, n)
798
t = cputime()
799
v = A.rank()
800
return cputime(t)
801
elif system == 'magma':
802
code = """
803
n := %s;
804
A := Random(MatrixAlgebra(GF(%s), n));
805
t := Cputime();
806
K := Rank(A);
807
s := Cputime(t);
808
"""%(n,p)
809
if verbose: print code
810
magma.eval(code)
811
return float(magma.eval('s'))
812
else:
813
raise ValueError, 'unknown system "%s"'%system
814
815
def det_GF(n=400, p=16411 , system='sage'):
816
"""
817
Dense determinant over GF(p).
818
Given an n x n matrix A over GF with random entries compute
819
det(A).
820
821
INPUT:
822
823
- ``n`` - matrix dimension (default: 300)
824
- ``p`` - prime number (default: ``16411``)
825
- ``system`` - either 'magma' or 'sage' (default: 'sage')
826
827
EXAMPLES::
828
829
sage: import sage.matrix.benchmark as b
830
sage: ts = b.det_GF(1000)
831
sage: tm = b.det_GF(1000, system='magma') # optional - magma
832
"""
833
if system == 'sage':
834
A = random_matrix(GF(p), n, n)
835
t = cputime()
836
d = A.determinant()
837
return cputime(t)
838
elif system == 'magma':
839
code = """
840
n := %s;
841
A := Random(MatrixAlgebra(GF(%s), n));
842
t := Cputime();
843
d := Determinant(A);
844
s := Cputime(t);
845
"""%(n,p)
846
if verbose: print code
847
magma.eval(code)
848
return float(magma.eval('s'))
849
else:
850
raise ValueError, 'unknown system "%s"'%system
851
852
853
#######################################################################
854
# Dense Benchmarks over QQ
855
#######################################################################
856
857
def hilbert_matrix(n):
858
"""
859
Returns the Hilbert matrix of size n over rationals.
860
861
EXAMPLES::
862
863
sage: import sage.matrix.benchmark as b
864
sage: b.hilbert_matrix(3)
865
[ 1 1/2 1/3]
866
[1/2 1/3 1/4]
867
[1/3 1/4 1/5]
868
"""
869
A = Matrix(QQ,n,n)
870
for i in range(A.nrows()):
871
for j in range(A.ncols()):
872
A[i,j] = QQ(1)/((i+1)+(j+1)-1)
873
return A
874
875
# Reduced row echelon form over QQ
876
877
def echelon_QQ(n=100, min=0, max=9, system='sage'):
878
"""
879
Given a n x (2*n) matrix over QQ with random integer entries
880
between min and max, compute the reduced row echelon form.
881
882
INPUT:
883
884
- ``n`` - matrix dimension (default: ``300``)
885
- ``min`` - minimal value for entries of matrix (default: ``-9``)
886
- ``max`` - maximal value for entries of matrix (default: ``9``)
887
- ``system`` - either 'sage' or 'magma' (default: 'sage')
888
889
EXAMPLES::
890
891
sage: import sage.matrix.benchmark as b
892
sage: ts = b.echelon_QQ(100)
893
sage: tm = b.echelon_QQ(100, system='magma') # optional - magma
894
"""
895
if system == 'sage':
896
A = random_matrix(ZZ, n, 2*n, x=min, y=max+1).change_ring(QQ)
897
t = cputime()
898
v = A.echelon_form()
899
return cputime(t)
900
elif system == 'magma':
901
code = """
902
n := %s;
903
A := RMatrixSpace(RationalField(), n, 2*n)![Random(%s,%s) : i in [1..n*2*n]];
904
t := Cputime();
905
K := EchelonForm(A);
906
s := Cputime(t);
907
"""%(n,min,max)
908
if verbose: print code
909
magma.eval(code)
910
return float(magma.eval('s'))
911
else:
912
raise ValueError, 'unknown system "%s"'%system
913
914
# Invert a matrix over QQ.
915
916
def inverse_QQ(n=100, min=0, max=9, system='sage'):
917
"""
918
Given a n x n matrix over QQ with random integer entries
919
between min and max, compute the reduced row echelon form.
920
921
INPUT:
922
923
- ``n`` - matrix dimension (default: ``300``)
924
- ``min`` - minimal value for entries of matrix (default: ``-9``)
925
- ``max`` - maximal value for entries of matrix (default: ``9``)
926
- ``system`` - either 'sage' or 'magma' (default: 'sage')
927
928
EXAMPLES::
929
930
sage: import sage.matrix.benchmark as b
931
sage: ts = b.inverse_QQ(100)
932
sage: tm = b.inverse_QQ(100, system='magma') # optional - magma
933
"""
934
if system == 'sage':
935
A = random_matrix(ZZ, n, n, x=min, y=max+1).change_ring(QQ)
936
t = cputime()
937
v = ~A
938
return cputime(t)
939
elif system == 'magma':
940
code = """
941
n := %s;
942
A := MatrixAlgebra(RationalField(), n)![Random(%s,%s) : i in [1..n*n]];
943
t := Cputime();
944
K := A^(-1);
945
s := Cputime(t);
946
"""%(n,min,max)
947
if verbose: print code
948
magma.eval(code)
949
return float(magma.eval('s'))
950
else:
951
raise ValueError, 'unknown system "%s"'%system
952
953
954
# Matrix multiplication over QQ
955
def matrix_multiply_QQ(n=100, bnd=2, system='sage', times=1):
956
"""
957
Given an n x n matrix A over QQ with random entries
958
whose numerators and denominators are bounded by bnd,
959
compute A * (A+1).
960
961
INPUT:
962
963
- ``n`` - matrix dimension (default: ``300``)
964
- ``bnd`` - numerator and denominator bound (default: ``bnd``)
965
- ``system`` - either 'sage' or 'magma' (default: 'sage')
966
- ``times`` - number of experiments (default: ``1``)
967
968
EXAMPLES::
969
970
sage: import sage.matrix.benchmark as b
971
sage: ts = b.matrix_multiply_QQ(100)
972
sage: tm = b.matrix_multiply_QQ(100, system='magma') # optional - magma
973
"""
974
if system == 'sage':
975
A = random_matrix(QQ, n, n, num_bound=bnd, den_bound=bnd)
976
B = A + 1
977
t = cputime()
978
for z in range(times):
979
v = A * B
980
return cputime(t)/times
981
elif system == 'magma':
982
A = magma(random_matrix(QQ, n, n, num_bound=bnd, den_bound=bnd))
983
code = """
984
n := %s;
985
A := %s;
986
B := A + 1;
987
t := Cputime();
988
for z in [1..%s] do
989
K := A * B;
990
end for;
991
s := Cputime(t);
992
"""%(n, A.name(), times)
993
if verbose: print code
994
magma.eval(code)
995
return float(magma.eval('s'))/times
996
else:
997
raise ValueError, 'unknown system "%s"'%system
998
999
1000
# Determinant of Hilbert matrix
1001
def det_hilbert_QQ(n=80, system='sage'):
1002
"""
1003
Runs the benchmark for calculating the determinant of the hilbert
1004
matrix over rationals of dimension n.
1005
1006
INPUT:
1007
1008
- ``n`` - matrix dimension (default: ``300``)
1009
- ``system`` - either 'sage' or 'magma' (default: 'sage')
1010
1011
EXAMPLES::
1012
1013
sage: import sage.matrix.benchmark as b
1014
sage: ts = b.det_hilbert_QQ(50)
1015
sage: tm = b.det_hilbert_QQ(50, system='magma') # optional - magma
1016
"""
1017
if system == 'sage':
1018
A = hilbert_matrix(n)
1019
t = cputime()
1020
d = A.determinant()
1021
return cputime(t)
1022
elif system == 'magma':
1023
code = """
1024
h := HilbertMatrix(%s);
1025
tinit := Cputime();
1026
d := Determinant(h);
1027
s := Cputime(tinit);
1028
delete h;
1029
"""%n
1030
if verbose: print code
1031
magma.eval(code)
1032
return float(magma.eval('s'))
1033
1034
# inverse of Hilbert matrix
1035
def invert_hilbert_QQ(n=40, system='sage'):
1036
"""
1037
Runs the benchmark for calculating the inverse of the hilbert
1038
matrix over rationals of dimension n.
1039
1040
INPUT:
1041
1042
- ``n`` - matrix dimension (default: ``300``)
1043
- ``system`` - either 'sage' or 'magma' (default: 'sage')
1044
1045
EXAMPLES::
1046
1047
sage: import sage.matrix.benchmark as b
1048
sage: ts = b.invert_hilbert_QQ(30)
1049
sage: tm = b.invert_hilbert_QQ(30, system='magma') # optional - magma
1050
"""
1051
if system == 'sage':
1052
A = hilbert_matrix(n)
1053
t = cputime()
1054
d = A**(-1)
1055
return cputime(t)
1056
elif system == 'magma':
1057
code = """
1058
h := HilbertMatrix(%s);
1059
tinit := Cputime();
1060
d := h^(-1);
1061
s := Cputime(tinit);
1062
delete h;
1063
"""%n
1064
if verbose: print code
1065
magma.eval(code)
1066
return float(magma.eval('s'))
1067
1068
def MatrixVector_QQ(n=1000,h=100,system='sage',times=1):
1069
"""
1070
Compute product of square ``n`` matrix by random vector with num and
1071
denom bounded by ``h`` the given number of ``times``.
1072
1073
INPUT:
1074
1075
- ``n`` - matrix dimension (default: ``300``)
1076
- ``h`` - numerator and denominator bound (default: ``bnd``)
1077
- ``system`` - either 'sage' or 'magma' (default: 'sage')
1078
- ``times`` - number of experiments (default: ``1``)
1079
1080
EXAMPLES::
1081
1082
sage: import sage.matrix.benchmark as b
1083
sage: ts = b.MatrixVector_QQ(500)
1084
sage: tm = b.MatrixVector_QQ(500, system='magma') # optional - magma
1085
"""
1086
if system=='sage':
1087
V=QQ**n
1088
v=V.random_element(h)
1089
M=random_matrix(QQ,n)
1090
t=cputime()
1091
for i in range(times):
1092
w=M*v
1093
return cputime(t)
1094
elif system == 'magma':
1095
code = """
1096
n:=%s;
1097
h:=%s;
1098
times:=%s;
1099
v:=VectorSpace(RationalField(),n)![Random(h)/(Random(h)+1) : i in [1..n]];
1100
M:=MatrixAlgebra(RationalField(),n)![Random(h)/(Random(h)+1) : i in [1..n^2]];
1101
t := Cputime();
1102
for z in [1..times] do
1103
W:=v*M;
1104
end for;
1105
s := Cputime(t);
1106
"""%(n,h,times)
1107
if verbose: print code
1108
magma.eval(code)
1109
return float(magma.eval('s'))
1110
else:
1111
raise ValueError, 'unknown system "%s"'%system
1112
1113
1114
#######################################################################
1115
# Dense Benchmarks over machine reals
1116
# Note that the precision in reals for MAGMA is base 10, while in
1117
# sage it is in base 2
1118
#######################################################################
1119
1120
# Real Nullspace
1121
1122
def nullspace_RR(n=300, min=0, max=10, system='sage'):
1123
"""
1124
Nullspace over RR:
1125
Given a n+1 x n matrix over RR with random entries
1126
between min and max, compute the nullspace.
1127
1128
INPUT:
1129
1130
- ``n`` - matrix dimension (default: ``300``)
1131
- ``min`` - minimal value for entries of matrix (default: ``0``)
1132
- ``max`` - maximal value for entries of matrix (default: ``10``)
1133
- ``system`` - either 'sage' or 'magma' (default: 'sage')
1134
1135
EXAMPLES::
1136
1137
sage: import sage.matrix.benchmark as b
1138
sage: ts = b.nullspace_RR(100)
1139
sage: tm = b.nullspace_RR(100, system='magma') # optional - magma
1140
"""
1141
if system == 'sage':
1142
A = random_matrix(ZZ, n+1, n, x=min, y=max+1).change_ring(RR)
1143
t = cputime()
1144
v = A.kernel()
1145
return cputime(t)
1146
elif system == 'magma':
1147
code = """
1148
n := %s;
1149
A := RMatrixSpace(RealField(16), n+1,n)![Random(%s,%s) : i in [1..n*(n+1)]];
1150
t := Cputime();
1151
K := Kernel(A);
1152
s := Cputime(t);
1153
"""%(n,min,max)
1154
if verbose: print code
1155
magma.eval(code)
1156
return float(magma.eval('s'))
1157
else:
1158
raise ValueError, 'unknown system "%s"'%system
1159
1160
1161
def nullspace_RDF(n=300, min=0, max=10, system='sage'):
1162
"""
1163
Nullspace over RDF:
1164
Given a n+1 x n matrix over RDF with random entries
1165
between min and max, compute the nullspace.
1166
1167
INPUT:
1168
1169
- ``n`` - matrix dimension (default: ``300``)
1170
- ``min`` - minimal value for entries of matrix (default: ``0``)
1171
- ``max`` - maximal value for entries of matrix (default: `10``)
1172
- ``system`` - either 'sage' or 'magma' (default: 'sage')
1173
1174
EXAMPLES::
1175
1176
sage: import sage.matrix.benchmark as b
1177
sage: ts = b.nullspace_RDF(100)
1178
sage: tm = b.nullspace_RDF(100, system='magma') # optional - magma
1179
"""
1180
if system == 'sage':
1181
A = random_matrix(ZZ, n+1, n, x=min, y=max+1).change_ring(RDF)
1182
t = cputime()
1183
v = A.kernel()
1184
return cputime(t)
1185
elif system == 'magma':
1186
code = """
1187
n := %s;
1188
A := RMatrixSpace(RealField(16), n+1,n)![Random(%s,%s) : i in [1..n*(n+1)]];
1189
t := Cputime();
1190
K := Kernel(A);
1191
s := Cputime(t);
1192
"""%(n,min,max)
1193
if verbose: print code
1194
magma.eval(code)
1195
return float(magma.eval('s'))
1196
else:
1197
raise ValueError, 'unknown system "%s"'%system
1198
1199
1200
1201