Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

r""" 

Linear feedback shift register (LFSR) sequence commands 

 

Stream ciphers have been used for a long time as a source of 

pseudo-random number generators. 

 

S. Golomb [Go1967]_ gives a list of three statistical properties a 

sequence of numbers `{\bf a}=\{a_n\}_{n=1}^\infty`, 

`a_n\in \{0,1\}`, should display to be considered 

"random". Define the autocorrelation of `{\bf a}` to be 

 

.. MATH:: 

 

C(k)=C(k,{\bf a})=\lim_{N\rightarrow \infty} {1\over N}\sum_{n=1}^N (-1)^{a_n+a_{n+k}}. 

 

 

In the case where `{\bf a}` is periodic with period 

`P` then this reduces to 

 

.. MATH:: 

 

C(k)={1\over P}\sum_{n=1}^P (-1)^{a_n+a_{n+k}}. 

 

 

Assume `{\bf a}` is periodic with period `P`. 

 

 

- balance: `|\sum_{n=1}^P(-1)^{a_n}|\leq 1`. 

 

- low autocorrelation: 

 

.. MATH:: 

 

C(k)= \left\{ \begin{array}{cc} 1,& k=0,\\ \epsilon, & k\not= 0. \end{array} \right. 

 

 

(For sequences satisfying these first two properties, it is known 

that `\epsilon=-1/P` must hold.) 

 

- proportional runs property: In each period, half the runs have 

length `1`, one-fourth have length `2`, etc. 

Moreover, there are as many runs of `1`'s as there are of 

`0`'s. 

 

 

A general feedback shift register is a map 

`f:{\bf F}_q^d\rightarrow {\bf F}_q^d` of the form 

 

.. MATH:: 

 

\begin{array}{c} f(x_0,...,x_{n-1})=(x_1,x_2,...,x_n),\\ x_n=C(x_0,...,x_{n-1}), \end{array} 

 

 

where `C:{\bf F}_q^d\rightarrow {\bf F}_q` is a given 

function. When `C` is of the form 

 

.. MATH:: 

 

C(x_0,...,x_{n-1})=a_0x_0+...+a_{n-1}x_{n-1}, 

 

 

for some given constants `a_i\in {\bf F}_q`, the map is 

called a linear feedback shift register (LFSR). 

 

Example of a LFSR Let 

 

.. MATH:: 

 

f(x)=a_{{0}}+a_{{1}}x+...+a_{{n}}{x}^n+..., 

 

 

 

.. MATH:: 

 

g(x)=b_{{0}}+b_{{1}}x+...+b_{{n}}{x}^n+..., 

 

 

be given polynomials in `{\bf F}_2[x]` and let 

 

.. MATH:: 

 

h(x)={f(x)\over g(x)}=c_0+c_1x+...+c_nx^n+... \ . 

 

 

We can compute a recursion formula which allows us to rapidly 

compute the coefficients of `h(x)` (take `f(x)=1`): 

 

.. MATH:: 

 

c_{n}=\sum_{i=1}^n {{-b_i\over b_0}c_{n-i}}. 

 

 

 

The coefficients of `h(x)` can, under certain conditions on 

`f(x)` and `g(x)`, be considered "random" from 

certain statistical points of view. 

 

Example: For instance, if 

 

.. MATH:: 

 

f(x)=1,\ \ \ \ g(x)=x^4+x+1, 

 

then 

 

.. MATH:: 

 

h(x)=1+x+x^2+x^3+x^5+x^7+x^8+...\ . 

 

The coefficients of `h` are 

 

.. MATH:: 

 

\begin{array}{c} 1, 1, 1, 0, 1, 0, 1, 1, 0, 0, 1, 0, 0, 0, 1, \\ 1, 1, 1, 0, 1, 0, 1, 1, 0, 0, 1, 0, 0, 0, 1, ...\ . \end{array} 

 

 

The sequence of `0,1`'s is periodic with period 

`P=2^4-1=15` and satisfies Golomb's three randomness 

conditions. However, this sequence of period 15 can be "cracked" 

(i.e., a procedure to reproduce `g(x)`) by knowing only 8 

terms! This is the function of the Berlekamp-Massey algorithm [Mas1969]_, 

implemented as ``berlekamp_massey.py``. 

 

AUTHORS: 

 

- Timothy Brock 

 

Created 11-24-2005 by wdj. Last updated 12-02-2005. 

""" 

 

########################################################################### 

# Copyright (C) 2006 Timothy Brock 

# and William Stein <wstein@gmail.com> 

# 

# Distributed under the terms of the GNU General Public License (GPL) 

# 

# http://www.gnu.org/licenses/ 

########################################################################### 

 

import copy 

 

from sage.structure.all import Sequence 

from sage.rings.all import Integer, PolynomialRing 

from sage.rings.finite_rings.finite_field_constructor import is_FiniteField 

 

 

def lfsr_sequence(key, fill, n): 

r""" 

This function creates an lfsr sequence. 

 

INPUT: 

 

 

- ``key`` - a list of finite field elements, 

[c_0,c_1,...,c_k]. 

 

- ``fill`` - the list of the initial terms of the lfsr 

sequence, [x_0,x_1,...,x_k]. 

 

- ``n`` - number of terms of the sequence that the 

function returns. 

 

 

OUTPUT: The lfsr sequence defined by 

`x_{n+1} = c_kx_n+...+c_0x_{n-k}`, for 

`n \leq k`. 

 

EXAMPLES:: 

 

sage: F = GF(2); l = F(1); o = F(0) 

sage: F = GF(2); S = LaurentSeriesRing(F,'x'); x = S.gen() 

sage: fill = [l,l,o,l]; key = [1,o,o,l]; n = 20 

sage: L = lfsr_sequence(key,fill,20); L 

[1, 1, 0, 1, 0, 1, 1, 0, 0, 1, 0, 0, 0, 1, 1, 1, 1, 0, 1, 0] 

sage: g = berlekamp_massey(L); g 

x^4 + x^3 + 1 

sage: (1)/(g.reverse()+O(x^20)) 

1 + x + x^2 + x^3 + x^5 + x^7 + x^8 + x^11 + x^15 + x^16 + x^17 + x^18 + O(x^20) 

sage: (1+x^2)/(g.reverse()+O(x^20)) 

1 + x + x^4 + x^8 + x^9 + x^10 + x^11 + x^13 + x^15 + x^16 + x^19 + O(x^20) 

sage: (1+x^2+x^3)/(g.reverse()+O(x^20)) 

1 + x + x^3 + x^5 + x^6 + x^9 + x^13 + x^14 + x^15 + x^16 + x^18 + O(x^20) 

sage: fill = [l,l,o,l]; key = [l,o,o,o]; n = 20 

sage: L = lfsr_sequence(key,fill,20); L 

[1, 1, 0, 1, 1, 1, 0, 1, 1, 1, 0, 1, 1, 1, 0, 1, 1, 1, 0, 1] 

sage: g = berlekamp_massey(L); g 

x^4 + 1 

sage: (1+x)/(g.reverse()+O(x^20)) 

1 + x + x^4 + x^5 + x^8 + x^9 + x^12 + x^13 + x^16 + x^17 + O(x^20) 

sage: (1+x+x^3)/(g.reverse()+O(x^20)) 

1 + x + x^3 + x^4 + x^5 + x^7 + x^8 + x^9 + x^11 + x^12 + x^13 + x^15 + x^16 + x^17 + x^19 + O(x^20) 

 

AUTHORS: 

 

- Timothy Brock (2005-11): with code modified from Python 

Cookbook, http://aspn.activestate.com/ASPN/Python/Cookbook/ 

""" 

if not isinstance(key, list): 

raise TypeError("key must be a list") 

key = Sequence(key) 

F = key.universe() 

if not is_FiniteField(F): 

raise TypeError("universe of sequence must be a finite field") 

 

s = fill 

k = len(fill) 

L = [] 

for i in range(n): 

s0 = copy.copy(s) 

L.append(s[0]) 

s = s[1:k] 

s.append(sum([key[i]*s0[i] for i in range(k)])) 

return L 

 

def lfsr_autocorrelation(L, p, k): 

""" 

INPUT: 

 

 

- ``L`` - is a periodic sequence of elements of ZZ or 

GF(2). L must have length = p 

 

- ``p`` - the period of L 

 

- ``k`` - k is an integer (0 k p) 

 

 

OUTPUT: autocorrelation sequence of L 

 

EXAMPLES:: 

 

sage: F = GF(2) 

sage: o = F(0) 

sage: l = F(1) 

sage: key = [l,o,o,l]; fill = [l,l,o,l]; n = 20 

sage: s = lfsr_sequence(key,fill,n) 

sage: lfsr_autocorrelation(s,15,7) 

4/15 

sage: lfsr_autocorrelation(s,int(15),7) 

4/15 

 

AUTHORS: 

 

- Timothy Brock (2006-04-17) 

""" 

if not isinstance(L, list): 

raise TypeError("L (=%s) must be a list"%L) 

p = Integer(p) 

_p = int(p) 

k = int(k) 

L0 = L[:_p] # slices makes a copy 

L0 = L0 + L0[:k] 

L1 = [int(L0[i])*int(L0[i + k])/p for i in range(_p)] 

return sum(L1) 

 

def lfsr_connection_polynomial(s): 

""" 

INPUT: 

 

 

- ``s`` - a sequence of elements of a finite field (F) 

of even length 

 

 

OUTPUT: 

 

 

- ``C(x)`` - the connection polynomial of the minimal 

LFSR. 

 

 

This implements the algorithm in section 3 of J. L. Massey's 

article [M]_. 

 

EXAMPLES:: 

 

sage: F = GF(2) 

sage: F 

Finite Field of size 2 

sage: o = F(0); l = F(1) 

sage: key = [l,o,o,l]; fill = [l,l,o,l]; n = 20 

sage: s = lfsr_sequence(key,fill,n); s 

[1, 1, 0, 1, 0, 1, 1, 0, 0, 1, 0, 0, 0, 1, 1, 1, 1, 0, 1, 0] 

sage: lfsr_connection_polynomial(s) 

x^4 + x + 1 

sage: berlekamp_massey(s) 

x^4 + x^3 + 1 

 

Notice that ``berlekamp_massey`` returns the reverse 

of the connection polynomial (and is potentially must faster than 

this implementation). 

 

AUTHORS: 

 

- Timothy Brock (2006-04-17) 

""" 

# Initialization: 

FF = s[0].base_ring() 

R = PolynomialRing(FF, "x") 

x = R.gen() 

C = R(1); B = R(1); m = 1; b = FF(1); L = 0; N = 0 

 

while N < len(s): 

if L > 0: 

r = min(L+1,C.degree()+1) 

d = s[N] + sum([(C.list())[i]*s[N-i] for i in range(1,r)]) 

if L == 0: 

d = s[N] 

if d == 0: 

m += 1 

N += 1 

if d > 0: 

if 2*L > N: 

C = C - d*b**(-1)*x**m*B 

m += 1 

N += 1 

else: 

T = C 

C = C - d*b**(-1)*x**m*B 

L = N + 1 - L 

m = 1 

b = d 

B = T 

N += 1 

return C