Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/crypto/lattice.py
4076 views
1
"""
2
Hard Lattice Generator
3
4
This module contains lattice related functions relevant in
5
cryptography.
6
7
Feel free to add more functionality.
8
9
AUTHORS:
10
* Richard Lindner <[email protected]>
11
* Michael Schneider <[email protected]>
12
"""
13
14
def gen_lattice(type='modular', n=4, m=8, q=11, seed=None, \
15
quotient=None, dual=False, ntl=False):
16
"""
17
This function generates different types of integral lattice bases
18
of row vectors relevant in cryptography.
19
20
Randomness can be set either with ``seed``, or by using
21
:func:`sage.misc.randstate.set_random_seed`.
22
23
INPUT:
24
25
* ``type`` - one of the following strings
26
* ``'modular'`` (default). A class of lattices for which
27
asymptotic worst-case to average-case connections hold. For
28
more refer to [A96]_.
29
* ``'random'`` - Special case of modular (n=1). A dense class
30
of lattice used for testing basis reduction algorithms
31
proposed by Goldstein and Mayer [GM02]_.
32
* ``'ideal'`` - Special case of modular. Allows for a more
33
compact representation proposed by [LM06]_.
34
* ``'cyclotomic'`` - Special case of ideal. Allows for
35
efficient processing proposed by [LM06]_.
36
* ``n`` - Determinant size, primal:`det(L) = q^n`, dual:`det(L) = q^{m-n}`.
37
For ideal lattices this is also the degree of the quotient polynomial.
38
* ``m`` - Lattice dimension, `L \subseteq Z^m`.
39
* ``q`` - Coefficent size, `q*Z^m \subseteq L`.
40
* ``seed`` - Randomness seed.
41
* ``quotient`` - For the type ideal, this determines the quotient
42
polynomial. Ignored for all other types.
43
* ``dual`` - Set this flag if you want a basis for `q*dual(L)`, for example
44
for Regev's LWE bases [R05]_.
45
* ``ntl`` - Set this flag if you want the lattice basis in NTL readable
46
format.
47
48
OUTPUT: ``B`` a unique size-reduced triangular (primal: lower_left,
49
dual: lower_right) basis of row vectors for the lattice in question.
50
51
EXAMPLES:
52
53
* Modular basis ::
54
55
sage: sage.crypto.gen_lattice(m=10, seed=42)
56
[11 0 0 0 0 0 0 0 0 0]
57
[ 0 11 0 0 0 0 0 0 0 0]
58
[ 0 0 11 0 0 0 0 0 0 0]
59
[ 0 0 0 11 0 0 0 0 0 0]
60
[ 2 4 3 5 1 0 0 0 0 0]
61
[ 1 -5 -4 2 0 1 0 0 0 0]
62
[-4 3 -1 1 0 0 1 0 0 0]
63
[-2 -3 -4 -1 0 0 0 1 0 0]
64
[-5 -5 3 3 0 0 0 0 1 0]
65
[-4 -3 2 -5 0 0 0 0 0 1]
66
67
* Random basis ::
68
69
sage: sage.crypto.gen_lattice(type='random', n=1, m=10, q=11^4, seed=42)
70
[14641 0 0 0 0 0 0 0 0 0]
71
[ 431 1 0 0 0 0 0 0 0 0]
72
[-4792 0 1 0 0 0 0 0 0 0]
73
[ 1015 0 0 1 0 0 0 0 0 0]
74
[-3086 0 0 0 1 0 0 0 0 0]
75
[-5378 0 0 0 0 1 0 0 0 0]
76
[ 4769 0 0 0 0 0 1 0 0 0]
77
[-1159 0 0 0 0 0 0 1 0 0]
78
[ 3082 0 0 0 0 0 0 0 1 0]
79
[-4580 0 0 0 0 0 0 0 0 1]
80
81
* Ideal bases with quotient x^n-1, m=2*n are NTRU bases ::
82
83
sage: sage.crypto.gen_lattice(type='ideal', seed=42, quotient=x^4-1)
84
[11 0 0 0 0 0 0 0]
85
[ 0 11 0 0 0 0 0 0]
86
[ 0 0 11 0 0 0 0 0]
87
[ 0 0 0 11 0 0 0 0]
88
[ 4 -2 -3 -3 1 0 0 0]
89
[-3 4 -2 -3 0 1 0 0]
90
[-3 -3 4 -2 0 0 1 0]
91
[-2 -3 -3 4 0 0 0 1]
92
93
* Cyclotomic bases with n=2^k are SWIFFT bases ::
94
95
sage: sage.crypto.gen_lattice(type='cyclotomic', seed=42)
96
[11 0 0 0 0 0 0 0]
97
[ 0 11 0 0 0 0 0 0]
98
[ 0 0 11 0 0 0 0 0]
99
[ 0 0 0 11 0 0 0 0]
100
[ 4 -2 -3 -3 1 0 0 0]
101
[ 3 4 -2 -3 0 1 0 0]
102
[ 3 3 4 -2 0 0 1 0]
103
[ 2 3 3 4 0 0 0 1]
104
105
* Dual modular bases are related to Regev's famous public-key
106
encryption [R05]_ ::
107
108
sage: sage.crypto.gen_lattice(type='modular', m=10, seed=42, dual=True)
109
[ 0 0 0 0 0 0 0 0 0 11]
110
[ 0 0 0 0 0 0 0 0 11 0]
111
[ 0 0 0 0 0 0 0 11 0 0]
112
[ 0 0 0 0 0 0 11 0 0 0]
113
[ 0 0 0 0 0 11 0 0 0 0]
114
[ 0 0 0 0 11 0 0 0 0 0]
115
[ 0 0 0 1 -5 -2 -1 1 -3 5]
116
[ 0 0 1 0 -3 4 1 4 -3 -2]
117
[ 0 1 0 0 -4 5 -3 3 5 3]
118
[ 1 0 0 0 -2 -1 4 2 5 4]
119
120
* Relation of primal and dual bases ::
121
122
sage: B_primal=sage.crypto.gen_lattice(m=10, q=11, seed=42)
123
sage: B_dual=sage.crypto.gen_lattice(m=10, q=11, seed=42, dual=True)
124
sage: B_dual_alt=transpose(11*B_primal.inverse()).change_ring(ZZ)
125
sage: B_dual_alt.hermite_form() == B_dual.hermite_form()
126
True
127
128
REFERENCES:
129
130
.. [A96] Miklos Ajtai.
131
Generating hard instances of lattice problems (extended abstract).
132
STOC, pp. 99--108, ACM, 1996.
133
134
.. [GM02] Daniel Goldstein and Andrew Mayer.
135
On the equidistribution of Hecke points.
136
Forum Mathematicum, 15:2, pp. 165--189, De Gruyter, 2003.
137
138
.. [LM06] Vadim Lyubashevsky and Daniele Micciancio.
139
Generalized compact knapsacks are collision resistant.
140
ICALP, pp. 144--155, Springer, 2006.
141
142
.. [R05] Oded Regev.
143
On lattices, learning with errors, random linear codes, and cryptography.
144
STOC, pp. 84--93, ACM, 2005.
145
"""
146
from sage.rings.finite_rings.integer_mod_ring \
147
import IntegerModRing
148
from sage.matrix.constructor import matrix, \
149
identity_matrix, block_matrix
150
from sage.matrix.matrix_space import MatrixSpace
151
from sage.rings.integer_ring import IntegerRing
152
if seed != None:
153
from sage.misc.randstate import set_random_seed
154
set_random_seed(seed)
155
156
if type == 'random':
157
if n != 1: raise ValueError('random bases require n = 1')
158
159
ZZ = IntegerRing()
160
ZZ_q = IntegerModRing(q)
161
A = identity_matrix(ZZ_q, n)
162
163
if type == 'random' or type == 'modular':
164
R = MatrixSpace(ZZ_q, m-n, n)
165
A = A.stack(R.random_element())
166
167
elif type == 'ideal':
168
if quotient == None: raise \
169
ValueError('ideal bases require a quotient polynomial')
170
x = quotient.default_variable()
171
if n != quotient.degree(x): raise \
172
ValueError('ideal bases require n = quotient.degree()')
173
R = ZZ_q[x].quotient(quotient, x)
174
for i in range(m//n):
175
A = A.stack(R.random_element().matrix())
176
177
elif type == 'cyclotomic':
178
from sage.rings.arith import euler_phi
179
from sage.misc.functional import cyclotomic_polynomial
180
181
# we assume that n+1 <= min( euler_phi^{-1}(n) ) <= 2*n
182
found = False
183
for k in range(2*n,n,-1):
184
if euler_phi(k) == n:
185
found = True
186
break
187
if not found: raise \
188
ValueError('cyclotomic bases require that n is an image of' + \
189
'Euler\'s totient function')
190
191
R = ZZ_q['x'].quotient(cyclotomic_polynomial(k, 'x'), 'x')
192
for i in range(m//n):
193
A = A.stack(R.random_element().matrix())
194
195
# switch from representatives 0,...,(q-1) to (1-q)/2,....,(q-1)/2
196
def minrep(a):
197
if abs(a-q) < abs(a): return a-q
198
else: return a
199
A_prime = A[n:m].lift().apply_map(minrep)
200
201
if not dual:
202
B = block_matrix([[ZZ(q), ZZ.zero()], [A_prime, ZZ.one()] ], \
203
subdivide=False)
204
else:
205
B = block_matrix([[ZZ.one(), -A_prime.transpose()], [ZZ.zero(), \
206
ZZ(q)]], subdivide=False)
207
for i in range(m//2): B.swap_rows(i,m-i-1)
208
209
if not ntl:
210
return B
211
else:
212
return B._ntl_()
213
214