Ticket #10848: trac_10848-hermitian-matrices.patch
File trac_10848-hermitian-matrices.patch, 6.0 KB (added by , 10 years ago) |
---|
-
sage/matrix/matrix0.pyx
# HG changeset patch # User Rob Beezer <beezer@ups.edu> # Date 1298574394 28800 # Node ID 0fb3e056a11a8ca32ac5b0325d0bdd04c4f10ded # Parent 3c67675a1df921b9a9fa0ef007fc54c1a1d1e8f6 10848: add is_hermitian for matrices diff -r 3c67675a1df9 -r 0fb3e056a11a sage/matrix/matrix0.pyx
a b 2896 2896 return False 2897 2897 return True 2898 2898 2899 def is_hermitian(self): 2900 r""" 2901 Returns ``True`` if the matrix is equal to its conjugate-transpose. 2902 2903 OUTPUT: 2904 2905 ``True`` if the matrix is square and equal to the tranpose 2906 with every entry conjugated, and ``False`` otherwise. 2907 2908 Note that if conjugation has no effect on elements of the base 2909 ring (such as for integers), then the :meth:`is_symmetric` 2910 method is equivalent and faster. 2911 2912 For numerical matrices a specialized routine available 2913 over ``RDF`` and ``CDF`` is a good choice. 2914 2915 EXAMPLES:: 2916 2917 sage: A = matrix(QQbar, [[ 1 + I, 1 - 6*I, -1 - I], 2918 ... [-3 - I, -4*I, -2], 2919 ... [-1 + I, -2 - 8*I, 2 + I]]) 2920 sage: A.is_hermitian() 2921 False 2922 sage: B = A*A.conjugate_transpose() 2923 sage: B.is_hermitian() 2924 True 2925 2926 Sage has several fields besides the entire complex numbers 2927 where conjugation is non-trivial. :: 2928 2929 sage: F.<b> = QuadraticField(-7) 2930 sage: C = matrix(F, [[-2*b - 3, 7*b - 6, -b + 3], 2931 ... [-2*b - 3, -3*b + 2, -2*b], 2932 ... [ b + 1, 0, -2]]) 2933 sage: C.is_hermitian() 2934 False 2935 sage: C = C*C.conjugate_transpose() 2936 sage: C.is_hermitian() 2937 True 2938 2939 Rectangular matrices are never Hermitian. :: 2940 2941 sage: A = matrix(QQbar, 3, 4) 2942 sage: A.is_hermitian() 2943 False 2944 """ 2945 if not self.is_square(): 2946 return False 2947 2948 cdef Py_ssize_t i,j 2949 2950 for i from 0 <= i < self._nrows: 2951 for j from 0 <= j < i: 2952 if self.get_unsafe(i,j) != self.get_unsafe(j,i).conjugate(): 2953 return False 2954 return True 2955 2899 2956 def is_skew_symmetric(self): 2900 2957 """ 2901 2958 Returns True if this is a skew symmetric matrix. -
sage/matrix/matrix_double_dense.pyx
diff -r 3c67675a1df9 -r 0fb3e056a11a sage/matrix/matrix_double_dense.pyx
a b 1327 1327 return False 1328 1328 b = True 1329 1329 for i from 0 < i < self._nrows: 1330 for j from 0 <= j < i: 1330 for j from 0 <= j < i: 1331 1331 if math.fabs(self.get_unsafe(i,j) - self.get_unsafe(j,i)) > tol: 1332 1332 b = False 1333 1333 break 1334 1334 self.cache(key, b) 1335 1335 return b 1336 1336 1337 def is_hermitian(self, tol = 1e-12): 1338 r""" 1339 Return ``True`` if the matrix is equal to its conjugate-transpose. 1340 1341 INPUT: 1342 1343 - ``tol`` - default: ``1e-12`` - the largest value of the 1344 absolute value of the difference between two matrix entries 1345 for which they will still be considered equal. 1346 1347 OUTPUT: 1348 1349 ``True`` if the matrix is square and equal to the tranpose 1350 with every entry conjugated, and ``False`` otherwise. 1351 1352 Note that if conjugation has no effect on elements of the base 1353 ring (such as for integers), then the :meth:`is_symmetric` 1354 method is equivalent and faster. 1355 1356 The tolerance parameter is used to allow for numerical values 1357 to be equal if there is a slight difference due to round-off 1358 and other imprecisions. 1359 1360 The result is cached, on a per-tolerance basis. 1361 1362 EXAMPLES:: 1363 1364 sage: A = matrix(CDF, [[ 1 + I, 1 - 6*I, -1 - I], 1365 ... [-3 - I, -4*I, -2], 1366 ... [-1 + I, -2 - 8*I, 2 + I]]) 1367 sage: A.is_hermitian() 1368 False 1369 sage: B = A*A.conjugate_transpose() 1370 sage: B.is_hermitian() 1371 True 1372 1373 We get a unitary matrix from the SVD routine and use this 1374 numerical matrix to create a matrix that should be Hermitian 1375 (indeed it should be the identity matrix), but with some 1376 imprecision. We use this to illustrate that if the tolerance 1377 is set too small, then we can be too strict about the equality 1378 of entries and achieve the wrong result. :: 1379 1380 sage: A = matrix(CDF, [[ 1 + I, 1 - 6*I, -1 - I], 1381 ... [-3 - I, -4*I, -2], 1382 ... [-1 + I, -2 - 8*I, 2 + I]]) 1383 sage: U, _, _ = A.SVD() 1384 sage: B=U*U.conjugate_transpose() 1385 sage: B.is_hermitian() 1386 True 1387 sage: B.is_hermitian(tol=5.0e-17) 1388 True 1389 sage: B.is_hermitian(tol=3.0e-17) 1390 False 1391 1392 Rectangular matrices are never Hermitian. :: 1393 1394 sage: A = matrix(CDF, 3, 4) 1395 sage: A.is_hermitian() 1396 False 1397 """ 1398 global numpy 1399 tol = float(tol) 1400 key = 'hermitian_%s'%tol 1401 b = self.fetch(key) 1402 if not b is None: 1403 return b 1404 if not self.is_square(): 1405 self.cache(key, False) 1406 return False 1407 if numpy is None: 1408 import numpy 1409 cdef Py_ssize_t i, j 1410 hermitian = True 1411 for i from 0 < i < self._nrows: 1412 for j from 0 <= j < i: 1413 if numpy.absolute(self.get_unsafe(i,j) - self.get_unsafe(j,i).conjugate()) > tol: 1414 hermitian = False 1415 break 1416 self.cache(key, hermitian) 1417 return hermitian 1337 1418 1338 1419 def cholesky(self): 1339 1420 r"""