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

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

r""" 

Orthogonal Linear Groups 

 

The general orthogonal group `GO(n,R)` consists of all `n \times n` 

matrices over the ring `R` preserving an `n`-ary positive definite 

quadratic form. In cases where there are multiple non-isomorphic 

quadratic forms, additional data needs to be specified to 

disambiguate. The special orthogonal group is the normal subgroup of 

matrices of determinant one. 

 

In characteristics different from 2, a quadratic form is equivalent to 

a bilinear symmetric form. Furthermore, over the real numbers a 

positive definite quadratic form is equivalent to the diagonal 

quadratic form, equivalent to the bilinear symmetric form defined by 

the identity matrix. Hence, the orthogonal group `GO(n,\RR)` is the 

group of orthogonal matrices in the usual sense. 

 

In the case of a finite field and if the degree `n` is even, then 

there are two inequivalent quadratic forms and a third parameter ``e`` 

must be specified to disambiguate these two possibilities. The index 

of `SO(e,d,q)` in `GO(e,d,q)` is `2` if `q` is odd, but `SO(e,d,q) = 

GO(e,d,q)` if `q` is even.) 

 

.. warning:: 

 

GAP and Sage use different notations: 

 

* GAP notation: The optional ``e`` comes first, that is, ``GO([e,] 

d, q)``, ``SO([e,] d, q)``. 

 

* Sage notation: The optional ``e`` comes last, the standard Python 

convention: ``GO(d, GF(q), e=0)``, ``SO( d, GF(q), e=0)``. 

 

EXAMPLES:: 

 

sage: GO(3,7) 

General Orthogonal Group of degree 3 over Finite Field of size 7 

 

sage: G = SO( 4, GF(7), 1); G 

Special Orthogonal Group of degree 4 and form parameter 1 over Finite Field of size 7 

sage: G.random_element() # random 

[4 3 5 2] 

[6 6 4 0] 

[0 4 6 0] 

[4 4 5 1] 

 

TESTS:: 

 

sage: G = GO(3, GF(5)) 

sage: latex(G) 

\text{GO}_{3}(\Bold{F}_{5}) 

sage: G = SO(3, GF(5)) 

sage: latex(G) 

\text{SO}_{3}(\Bold{F}_{5}) 

sage: G = SO(4, GF(5), 1) 

sage: latex(G) 

\text{SO}_{4}(\Bold{F}_{5}, +) 

 

AUTHORS: 

 

- David Joyner (2006-03): initial version 

 

- David Joyner (2006-05): added examples, _latex_, __str__, gens, 

as_matrix_group 

 

- William Stein (2006-12-09): rewrite 

 

- Volker Braun (2013-1) port to new Parent, libGAP, extreme refactoring. 

""" 

 

#***************************************************************************** 

# Copyright (C) 2006 David Joyner and William Stein 

# Copyright (C) 2013 Volker Braun <vbraun.name@gmail.com> 

# 

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

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

#***************************************************************************** 

 

from sage.rings.all import ZZ 

from sage.rings.finite_rings.finite_field_base import is_FiniteField 

from sage.misc.latex import latex 

from sage.misc.cachefunc import cached_method 

from sage.groups.matrix_gps.named_group import ( 

normalize_args_vectorspace, NamedMatrixGroup_generic, NamedMatrixGroup_gap ) 

 

def normalize_args_e(degree, ring, e): 

""" 

Normalize the arguments that relate the choice of quadratic form 

for special orthogonal groups over finite fields. 

 

INPUT: 

 

- ``degree`` -- integer. The degree of the affine group, that is, 

the dimension of the affine space the group is acting on. 

 

- ``ring`` -- a ring. The base ring of the affine space. 

 

- ``e`` -- integer, one of `+1`, `0`, `-1`. Only relevant for 

finite fields and if the degree is even. A parameter that 

distinguishes inequivalent invariant forms. 

 

OUTPUT: 

 

The integer ``e`` with values required by GAP. 

 

TESTS:: 

 

sage: from sage.groups.matrix_gps.orthogonal import normalize_args_e 

sage: normalize_args_e(2, GF(3), +1) 

1 

sage: normalize_args_e(3, GF(3), 0) 

0 

sage: normalize_args_e(3, GF(3), +1) 

0 

sage: normalize_args_e(2, GF(3), 0) 

Traceback (most recent call last): 

... 

ValueError: must have e=-1 or e=1 for even degree 

""" 

if is_FiniteField(ring) and degree%2 == 0: 

if e not in (-1, +1): 

raise ValueError('must have e=-1 or e=1 for even degree') 

else: 

e = 0 

return ZZ(e) 

 

 

 

 

 

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

# General Orthogonal Group 

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

 

def GO(n, R, e=0, var='a'): 

""" 

Return the general orthogonal group. 

 

The general orthogonal group `GO(n,R)` consists of all `n \times n` 

matrices over the ring `R` preserving an `n`-ary positive definite 

quadratic form. In cases where there are multiple non-isomorphic 

quadratic forms, additional data needs to be specified to 

disambiguate. 

 

In the case of a finite field and if the degree `n` is even, then 

there are two inequivalent quadratic forms and a third parameter 

``e`` must be specified to disambiguate these two possibilities. 

 

.. note:: 

 

This group is also available via ``groups.matrix.GO()``. 

 

INPUT: 

 

- ``n`` -- integer. The degree. 

 

- ``R`` -- ring or an integer. If an integer is specified, the 

corresponding finite field is used. 

 

- ``e`` -- ``+1`` or ``-1``, and ignored by default. Only relevant 

for finite fields and if the degree is even. A parameter that 

distinguishes inequivalent invariant forms. 

 

OUTPUT: 

 

The general orthogonal group of given degree, base ring, and 

choice of invariant form. 

 

EXAMPLES:: 

 

sage: GO( 3, GF(7)) 

General Orthogonal Group of degree 3 over Finite Field of size 7 

sage: GO( 3, GF(7)).order() 

672 

sage: GO( 3, GF(7)).gens() 

( 

[3 0 0] [0 1 0] 

[0 5 0] [1 6 6] 

[0 0 1], [0 2 1] 

) 

 

TESTS:: 

 

sage: groups.matrix.GO(2, 3, e=-1) 

General Orthogonal Group of degree 2 and form parameter -1 over Finite Field of size 3 

""" 

degree, ring = normalize_args_vectorspace(n, R, var=var) 

e = normalize_args_e(degree, ring, e) 

if e == 0: 

name = 'General Orthogonal Group of degree {0} over {1}'.format(degree, ring) 

ltx = r'\text{{GO}}_{{{0}}}({1})'.format(degree, latex(ring)) 

else: 

name = 'General Orthogonal Group of degree' + \ 

' {0} and form parameter {1} over {2}'.format(degree, e, ring) 

ltx = r'\text{{GO}}_{{{0}}}({1}, {2})'.format(degree, latex(ring), '+' if e == 1 else '-') 

if is_FiniteField(ring): 

cmd = 'GO({0}, {1}, {2})'.format(e, degree, ring.characteristic()) 

return OrthogonalMatrixGroup_gap(degree, ring, False, name, ltx, cmd) 

else: 

return OrthogonalMatrixGroup_generic(degree, ring, False, name, ltx) 

 

 

 

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

# Special Orthogonal Group 

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

 

def SO(n, R, e=None, var='a'): 

""" 

Return the special orthogonal group. 

 

The special orthogonal group `GO(n,R)` consists of all `n \times n` 

matrices with determinant one over the ring `R` preserving an 

`n`-ary positive definite quadratic form. In cases where there are 

multiple non-isomorphic quadratic forms, additional data needs to 

be specified to disambiguate. 

 

.. note:: 

 

This group is also available via ``groups.matrix.SO()``. 

 

INPUT: 

 

- ``n`` -- integer. The degree. 

 

- ``R`` -- ring or an integer. If an integer is specified, the 

corresponding finite field is used. 

 

- ``e`` -- ``+1`` or ``-1``, and ignored by default. Only relevant 

for finite fields and if the degree is even. A parameter that 

distinguishes inequivalent invariant forms. 

 

OUTPUT: 

 

The special orthogonal group of given degree, base ring, and choice of 

invariant form. 

 

EXAMPLES:: 

 

sage: G = SO(3,GF(5)) 

sage: G 

Special Orthogonal Group of degree 3 over Finite Field of size 5 

 

sage: G = SO(3,GF(5)) 

sage: G.gens() 

( 

[2 0 0] [3 2 3] [1 4 4] 

[0 3 0] [0 2 0] [4 0 0] 

[0 0 1], [0 3 1], [2 0 4] 

) 

sage: G = SO(3,GF(5)) 

sage: G.as_matrix_group() 

Matrix group over Finite Field of size 5 with 3 generators ( 

[2 0 0] [3 2 3] [1 4 4] 

[0 3 0] [0 2 0] [4 0 0] 

[0 0 1], [0 3 1], [2 0 4] 

) 

 

TESTS:: 

 

sage: groups.matrix.SO(2, 3, e=1) 

Special Orthogonal Group of degree 2 and form parameter 1 over Finite Field of size 3 

""" 

degree, ring = normalize_args_vectorspace(n, R, var=var) 

e = normalize_args_e(degree, ring, e) 

if e == 0: 

name = 'Special Orthogonal Group of degree {0} over {1}'.format(degree, ring) 

ltx = r'\text{{SO}}_{{{0}}}({1})'.format(degree, latex(ring)) 

else: 

name = 'Special Orthogonal Group of degree' + \ 

' {0} and form parameter {1} over {2}'.format(degree, e, ring) 

ltx = r'\text{{SO}}_{{{0}}}({1}, {2})'.format(degree, latex(ring), '+' if e == 1 else '-') 

if is_FiniteField(ring): 

cmd = 'SO({0}, {1}, {2})'.format(e, degree, ring.characteristic()) 

return OrthogonalMatrixGroup_gap(degree, ring, True, name, ltx, cmd) 

else: 

return OrthogonalMatrixGroup_generic(degree, ring, True, name, ltx) 

 

 

 

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

# Orthogonal Group class 

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

 

class OrthogonalMatrixGroup_generic(NamedMatrixGroup_generic): 

 

@cached_method 

def invariant_bilinear_form(self): 

""" 

Return the symmetric bilinear form preserved by the orthogonal 

group. 

 

OUTPUT: 

 

A matrix. 

 

EXAMPLES:: 

 

sage: GO(2,3,+1).invariant_bilinear_form() 

[0 1] 

[1 0] 

sage: GO(2,3,-1).invariant_bilinear_form() 

[2 1] 

[1 1] 

""" 

from sage.matrix.constructor import identity_matrix 

m = identity_matrix(self.base_ring(), self.degree()) 

m.set_immutable() 

return m 

 

@cached_method 

def invariant_quadratic_form(self): 

""" 

Return the quadratic form preserved by the orthogonal group. 

 

OUTPUT: 

 

A matrix. 

 

EXAMPLES:: 

 

sage: GO(2,3,+1).invariant_quadratic_form() 

[0 1] 

[0 0] 

sage: GO(2,3,-1).invariant_quadratic_form() 

[1 1] 

[0 2] 

""" 

from sage.matrix.constructor import identity_matrix 

m = identity_matrix(self.base_ring(), self.degree()) 

m.set_immutable() 

return m 

 

def _check_matrix(self, x, *args): 

"""a 

Check whether the matrix ``x`` is symplectic. 

 

See :meth:`~sage.groups.matrix_gps.matrix_group._check_matrix` 

for details. 

 

EXAMPLES:: 

 

sage: G = GO(4, GF(5), +1) 

sage: G._check_matrix(G.an_element().matrix()) 

""" 

if self._special and x.determinant() != 1: 

raise TypeError('matrix must have determinant one') 

F = self.invariant_bilinear_form() 

if x * F * x.transpose() != F: 

raise TypeError('matrix must be orthogonal with respect to the invariant form') 

# TODO: check that quadratic form is preserved in characteristic two 

 

class OrthogonalMatrixGroup_gap(OrthogonalMatrixGroup_generic, NamedMatrixGroup_gap): 

 

@cached_method 

def invariant_bilinear_form(self): 

""" 

Return the symmetric bilinear form preserved by the orthogonal 

group. 

 

OUTPUT: 

 

A matrix `M` such that, for every group element g, the 

identity `g m g^T = m` holds. In characteristic different from 

two, this uniquely determines the orthogonal group. 

 

EXAMPLES:: 

 

sage: G = GO(4, GF(7), -1) 

sage: G.invariant_bilinear_form() 

[0 1 0 0] 

[1 0 0 0] 

[0 0 2 0] 

[0 0 0 2] 

 

sage: G = GO(4, GF(7), +1) 

sage: G.invariant_bilinear_form() 

[0 1 0 0] 

[1 0 0 0] 

[0 0 6 0] 

[0 0 0 2] 

 

sage: G = GO(4, QQ) 

sage: G.invariant_bilinear_form() 

[1 0 0 0] 

[0 1 0 0] 

[0 0 1 0] 

[0 0 0 1] 

 

sage: G = SO(4, GF(7), -1) 

sage: G.invariant_bilinear_form() 

[0 1 0 0] 

[1 0 0 0] 

[0 0 2 0] 

[0 0 0 2] 

""" 

m = self.gap().InvariantBilinearForm()['matrix'].matrix() 

m.set_immutable() 

return m 

 

@cached_method 

def invariant_quadratic_form(self): 

""" 

Return the quadratic form preserved by the orthogonal group. 

 

OUTPUT: 

 

The matrix `Q` defining "orthogonal" as follows. The matrix 

determines a quadratic form `q` on the natural vector space 

`V`, on which `G` acts, by `q(v) = v Q v^t`. A matrix `M` is 

an element of the orthogonal group if `q(v) = q(v M)` for all 

`v \in V`. 

 

EXAMPLES:: 

 

sage: G = GO(4, GF(7), -1) 

sage: G.invariant_quadratic_form() 

[0 1 0 0] 

[0 0 0 0] 

[0 0 1 0] 

[0 0 0 1] 

 

sage: G = GO(4, GF(7), +1) 

sage: G.invariant_quadratic_form() 

[0 1 0 0] 

[0 0 0 0] 

[0 0 3 0] 

[0 0 0 1] 

 

sage: G = GO(4, QQ) 

sage: G.invariant_quadratic_form() 

[1 0 0 0] 

[0 1 0 0] 

[0 0 1 0] 

[0 0 0 1] 

 

sage: G = SO(4, GF(7), -1) 

sage: G.invariant_quadratic_form() 

[0 1 0 0] 

[0 0 0 0] 

[0 0 1 0] 

[0 0 0 1] 

""" 

m = self.gap().InvariantQuadraticForm()['matrix'].matrix() 

m.set_immutable() 

return m