def iv_P_norm_expm(P_sqrt_T, M1, A, M2, tau): """ Bound on P-ellipsoid norm of (M1 (expm(A*t) - I) M2) for |t| < tau using the theorem in arXiv:1911.02537, section "Norm bounding of summands" @param P_sqrt_T: see iv_P_norm() """ P_sqrt_T = iv.matrix(P_sqrt_T) M1 = iv.matrix(M1) A = iv.matrix(A) M2 = iv.matrix(M2) # coerce tau to maximum tau = abs(iv.mpf(tau)).b # P-ellipsoid norms M1_p = iv_P_norm(M=M1, P_sqrt_T=P_sqrt_T) M2_p = iv_P_norm(M=M2, P_sqrt_T=P_sqrt_T) A_p = iv_P_norm(M=A, P_sqrt_T=P_sqrt_T) # A_pow[i] = A ** i A_pow = _iv_matrix_powers(A) # Work around bug in mpmath, see comment in iv_P_norm() zero = iv.matrix(mp.zeros(len(A))) M1 = zero + M1 # terms from [arXiv:1911.02537] M1_Ai_M2_p = lambda i: iv_P_norm(M=M1 @ A_pow[i] @ M2, P_sqrt_T=P_sqrt_T) gamma = lambda i: 1 / math.factorial(i) * (M1_Ai_M2_p(i) - M1_p * A_p ** i * M2_p) max_norm = sum([gamma(i) * (tau ** i) for i in range(1, IV_NORM_EVAL_ORDER + 1)]) + M1_p * M2_p * (iv.exp(A_p * tau) - 1) # the lower bound is always 0 (for t=0) return mp.mpi([0, max_norm.b])
def iv_P_norm(M, P_sqrt_T): """ interval bound on P-ellipsoid norm of M, defined as max_{x in R^n} sqrt(((M x).T P (M x)) / (x.T P x)) with P_sqrt_T.T @ P_sqrt_T = P, where x.T P x typically is a Lyapunov function """ P_sqrt_T = iv.matrix(P_sqrt_T) M = iv.matrix(M) # P_norm(M) = spectral norm (maximum singular value) of P_sqrt_T A P_sqrt_T**(-1) return iv_spectral_norm(iv.matrix(P_sqrt_T) @ M @ iv_matrix_inv(iv.matrix(P_sqrt_T)))
def example_matrices(self, include_singular=True, random=100): examples = [ iv.matrix([[1, 2], [4, 5]]) + iv.mpf([-1, +1]) * mp.mpf(1e-10), iv.eye(2) * 0.5, iv.diag([1, 1e-10]), iv.eye(2) * 1e-302 ] for n in [1, 4, 17]: for i in range(random // n): examples.append(iv.matrix(mp.randmatrix(n))) if include_singular: examples.append(iv.matrix(mp.zeros(4))) examples.append(iv.diag([1, 1e-200])) return examples
def iv_spectral_norm(M): """ Good interval bound of spectral norm of a (interval) matrix Theorem 3.2 from Siegfried M. Rump. “Verified bounds for singular values, in particular for the spectral norm of a matrix and its inverse”. In: BIT Numerical Mathematics 51.2 (Nov. 2010), pp. 367–384. DOI : 10.1007/s10543-010-0294-0. """ M = iv.matrix(M) # imprecise SVD of M (no requirement on precision) # _, _, V_T = mp.svd(iv_matrix_mid_to_mp(M)) # <-- configurable precision _, _, V_T = scipy.linalg.svd(iv_matrix_mid_to_numpy_ndarray(M)) # <-- faster # in [Rump 2010], SVD is defined as M = U @ S @ V.T, # in mpmath it is M = U @ S @ V (not transposed) for U,S,V=svd(M) # in scipy it is effectively the same as in mpmath, M = U @ S @ Vh for U,S,Vh = svd(M) V = numpy_ndarray_to_mp_matrix(V_T).T # now, everything is named as in [Rump2010], except that here, A is called M. # all following computations are interval bounds V = iv.matrix(V) B = M @ V BTB = B.T @ B # split BTB such that DE = diagonal D + rest E D = BTB @ 0 E = BTB @ 0 for i in range(BTB.rows): for j in range(BTB.cols): if i == j: D[i,j] = BTB[i,j] else: E[i,j] = BTB[i,j] # upper bound of spectral norm of I - V.T @ V alpha = iv_spectral_norm_rough(iv.eye(len(M)) - V.T @ V) # upper bound of spectral norm of E epsilon = iv_spectral_norm_rough(E) # maximum of D[i,i] (which are always >= 0) d_max = iv.norm(D, mp.inf) if alpha.b >= 1: # this shouldn't happen - even an imprecise SVD will roughly have V.T @ V = I. raise scipy.linalg.LinAlgError("Something's numerically wrong - the singular vectors are far from orthonormal") # should this ever happen in reality, a valid return value would be: # return iv.mpf([0, mp.inf]) try: lower_bound = iv.sqrt((d_max - epsilon) / (1 + alpha)).a except mp.libmp.libmpf.ComplexResult: lower_bound=0; # note that d_max, epsilon,alpha are intervals, so everything in the following computation is interval arithmetic return iv.mpf([lower_bound, iv.sqrt((d_max + epsilon) / (1 - alpha)).b])
def test_P_synthesis(self): """ Combined test: For A with eigenvalues inside the unit disk, generate P_sqrt_T such that P_norm(A, P_sqrt_T) < 1. """ for c in [0.001234, 1, -2, 42.123]: eigenvalue = 0.99 # rotation matrix with eigenvalue magnitude 0.99, # transformed with factor c (c=1: invariant set is a circle, otherwise: invariant set is elliptical) A = iv.matrix([[0., -c * eigenvalue], [eigenvalue / c, 0]]) # Note that A must be well-conditioned, otherwise this test will fail # compute eigenvalues of A eigv_A, _ = mp.eig(iv_matrix_mid_as_mp(A)) self.assertAlmostEqual(eigenvalue, max([abs(i) for i in eigv_A]), 3) P_sqrt_T = approx_P_sqrt_T(A) P_norm_iv = iv_P_norm(A, P_sqrt_T) self.check_P_norm(A, P_sqrt_T) # P_norm_iv must be at least 0.99, but should not be much larger than that self.assertLess(P_norm_iv, 1) self.assertGreater(P_norm_iv.b, 0.99) self.check_P_norm_expm(P_sqrt_T, M1=mp.randmatrix(2), A=A, M2=mp.randmatrix(2), tau=0.01)
def iv_matrix_inv(M): """ interval matrix inverse, with caching NOTE: If you change mpmath's resolution after the first call to this function, you may get cached old output with the old resolution. """ assert isinstance(M, (mp.matrix, iv.matrix)) M = iv.matrix(M) return M**(-1)
def _iv_matrix_powers(A): """ return the first IV_NORM_EVAL_ORDER+1 powers of A: [I, A, A**2, ..., A**(IV_NORM_EVAL_ORDER)] """ assert isinstance(A, iv.matrix) A = iv.matrix(A) + mp.mpi(0,0) # workaround bug A_pow = [iv.eye(len(A))] for i in range(IV_NORM_EVAL_ORDER+1): A_pow.append(A @ A_pow[-1]) return A_pow
def example_C(): """ DigitalControlLoop and P_sqrt_T value for example C. """ s = examples.example_C_quadrotor_attitude_one_axis() # P_sqrt_T = analyze(s, datatype=np)['P_sqrt_T'] P_sqrt_T = iv.matrix( [[ '0.0019244615174829', '0.010762303718393', '-0.000842841208611818', '-0.000546543030195919', '0.337150563056583' ], [ '0.0', '0.0145387640958695', '0.000181643285522379', '0.000217863576839351', '0.414793996245936' ], [ '0.0', '0.0', '0.443015131922611', '-0.443014979405948', '0.000266278043438505' ], ['0.0', '0.0', '0.0', '0.000392820699869412', '0.432006259021524'], ['0.0', '0.0', '0.0', '0.0', '0.662735066069117']]) return (s, P_sqrt_T)
def convert(M): if isinstance(M, mpmath.matrix): return iv.matrix(M) else: return iv.matrix(M.tolist())
def test_matrix_abs_max(self): A = iv.matrix([[1, 2], [4, 5]]) + iv.mpf([-1, +1]) * mp.mpf(1e-10) assert mp.matrix([[1, 2], [4, 5]]) + mp.mpf(1e-10) == matrix_abs_max(A) assert mp.matrix([[1, 2], [4, 5] ]) + mp.mpf(1e-10) == matrix_abs_max(-A)