brians | 343bc11 | 2013-02-10 01:53:46 +0000 | [diff] [blame] | 1 | #!/usr/bin/python |
| 2 | |
| 3 | """ |
| 4 | Control loop pole placement library. |
| 5 | |
| 6 | This library will grow to support many different pole placement methods. |
| 7 | Currently it only supports direct pole placement. |
| 8 | """ |
| 9 | |
| 10 | __author__ = 'Austin Schuh (austin.linux@gmail.com)' |
| 11 | |
| 12 | import numpy |
| 13 | import slycot |
Austin Schuh | c976f49 | 2015-02-22 21:28:18 -0800 | [diff] [blame] | 14 | import scipy.signal.cont2discrete |
Austin Schuh | c9177b5 | 2015-11-28 13:18:31 -0800 | [diff] [blame] | 15 | import glog |
brians | 343bc11 | 2013-02-10 01:53:46 +0000 | [diff] [blame] | 16 | |
| 17 | class Error (Exception): |
| 18 | """Base class for all control loop exceptions.""" |
| 19 | |
| 20 | |
| 21 | class PolePlacementError(Error): |
| 22 | """Exception raised when pole placement fails.""" |
| 23 | |
| 24 | |
| 25 | # TODO(aschuh): dplace should take a control system object. |
| 26 | # There should also exist a function to manipulate laplace expressions, and |
| 27 | # something to plot bode plots and all that. |
| 28 | def dplace(A, B, poles, alpha=1e-6): |
| 29 | """Set the poles of (A - BF) to poles. |
| 30 | |
| 31 | Args: |
| 32 | A: numpy.matrix(n x n), The A matrix. |
| 33 | B: numpy.matrix(n x m), The B matrix. |
| 34 | poles: array(imaginary numbers), The poles to use. Complex conjugates poles |
| 35 | must be in pairs. |
| 36 | |
| 37 | Raises: |
| 38 | ValueError: Arguments were the wrong shape or there were too many poles. |
| 39 | PolePlacementError: Pole placement failed. |
| 40 | |
| 41 | Returns: |
| 42 | numpy.matrix(m x n), K |
| 43 | """ |
| 44 | # See http://www.icm.tu-bs.de/NICONET/doc/SB01BD.html for a description of the |
| 45 | # fortran code that this is cleaning up the interface to. |
| 46 | n = A.shape[0] |
| 47 | if A.shape[1] != n: |
| 48 | raise ValueError("A must be square") |
| 49 | if B.shape[0] != n: |
| 50 | raise ValueError("B must have the same number of states as A.") |
| 51 | m = B.shape[1] |
| 52 | |
| 53 | num_poles = len(poles) |
| 54 | if num_poles > n: |
| 55 | raise ValueError("Trying to place more poles than states.") |
| 56 | |
| 57 | out = slycot.sb01bd(n=n, |
| 58 | m=m, |
| 59 | np=num_poles, |
| 60 | alpha=alpha, |
| 61 | A=A, |
| 62 | B=B, |
| 63 | w=numpy.array(poles), |
| 64 | dico='D') |
| 65 | |
| 66 | A_z = numpy.matrix(out[0]) |
| 67 | num_too_small_eigenvalues = out[2] |
| 68 | num_assigned_eigenvalues = out[3] |
| 69 | num_uncontrollable_eigenvalues = out[4] |
| 70 | K = numpy.matrix(-out[5]) |
| 71 | Z = numpy.matrix(out[6]) |
| 72 | |
| 73 | if num_too_small_eigenvalues != 0: |
| 74 | raise PolePlacementError("Number of eigenvalues that are too small " |
| 75 | "and are therefore unmodified is %d." % |
| 76 | num_too_small_eigenvalues) |
| 77 | if num_assigned_eigenvalues != num_poles: |
| 78 | raise PolePlacementError("Did not place all the eigenvalues that were " |
| 79 | "requested. Only placed %d eigenvalues." % |
| 80 | num_assigned_eigenvalues) |
| 81 | if num_uncontrollable_eigenvalues != 0: |
| 82 | raise PolePlacementError("Found %d uncontrollable eigenvlaues." % |
| 83 | num_uncontrollable_eigenvalues) |
| 84 | |
| 85 | return K |
Austin Schuh | c8ca244 | 2013-02-23 12:29:33 -0800 | [diff] [blame] | 86 | |
| 87 | |
| 88 | def c2d(A, B, dt): |
| 89 | """Converts from continuous time state space representation to discrete time. |
Austin Schuh | c8ca244 | 2013-02-23 12:29:33 -0800 | [diff] [blame] | 90 | Returns (A, B). C and D are unchanged.""" |
Austin Schuh | c8ca244 | 2013-02-23 12:29:33 -0800 | [diff] [blame] | 91 | |
Austin Schuh | c6d7a0f | 2017-02-04 22:15:20 -0800 | [diff] [blame] | 92 | ans_a, ans_b, _, _, _ = scipy.signal.cont2discrete( |
| 93 | (numpy.array(A), numpy.array(B), None, None), dt) |
Austin Schuh | c976f49 | 2015-02-22 21:28:18 -0800 | [diff] [blame] | 94 | return numpy.matrix(ans_a), numpy.matrix(ans_b) |
Austin Schuh | 7ec34fd | 2014-02-15 22:27:46 -0800 | [diff] [blame] | 95 | |
| 96 | def ctrb(A, B): |
Brian Silverman | e18cf50 | 2015-11-28 23:04:43 +0000 | [diff] [blame] | 97 | """Returns the controllability matrix. |
Austin Schuh | 7ec34fd | 2014-02-15 22:27:46 -0800 | [diff] [blame] | 98 | |
Austin Schuh | c9177b5 | 2015-11-28 13:18:31 -0800 | [diff] [blame] | 99 | This matrix must have full rank for all the states to be controllable. |
Austin Schuh | 7ec34fd | 2014-02-15 22:27:46 -0800 | [diff] [blame] | 100 | """ |
| 101 | n = A.shape[0] |
| 102 | output = B |
| 103 | intermediate = B |
| 104 | for i in xrange(0, n): |
| 105 | intermediate = A * intermediate |
| 106 | output = numpy.concatenate((output, intermediate), axis=1) |
| 107 | |
| 108 | return output |
| 109 | |
| 110 | def dlqr(A, B, Q, R): |
| 111 | """Solves for the optimal lqr controller. |
| 112 | |
| 113 | x(n+1) = A * x(n) + B * u(n) |
| 114 | J = sum(0, inf, x.T * Q * x + u.T * R * u) |
| 115 | """ |
| 116 | |
| 117 | # P = (A.T * P * A) - (A.T * P * B * numpy.linalg.inv(R + B.T * P *B) * (A.T * P.T * B).T + Q |
| 118 | |
Austin Schuh | 1a38796 | 2015-01-31 16:36:20 -0800 | [diff] [blame] | 119 | P, rcond, w, S, T = slycot.sb02od( |
| 120 | n=A.shape[0], m=B.shape[1], A=A, B=B, Q=Q, R=R, dico='D') |
Austin Schuh | 7ec34fd | 2014-02-15 22:27:46 -0800 | [diff] [blame] | 121 | |
| 122 | F = numpy.linalg.inv(R + B.T * P *B) * B.T * P * A |
| 123 | return F |
Austin Schuh | e4a14f2 | 2015-03-01 00:12:29 -0800 | [diff] [blame] | 124 | |
| 125 | def kalman(A, B, C, Q, R): |
| 126 | """Solves for the steady state kalman gain and covariance matricies. |
| 127 | |
| 128 | Args: |
| 129 | A, B, C: SS matricies. |
| 130 | Q: The model uncertantity |
| 131 | R: The measurement uncertainty |
| 132 | |
| 133 | Returns: |
| 134 | KalmanGain, Covariance. |
| 135 | """ |
Austin Schuh | 572ff40 | 2015-11-08 12:17:50 -0800 | [diff] [blame] | 136 | I = numpy.matrix(numpy.eye(Q.shape[0])) |
| 137 | Z = numpy.matrix(numpy.zeros(Q.shape[0])) |
Austin Schuh | c9177b5 | 2015-11-28 13:18:31 -0800 | [diff] [blame] | 138 | n = A.shape[0] |
| 139 | m = C.shape[0] |
| 140 | |
| 141 | controllability_rank = numpy.linalg.matrix_rank(ctrb(A.T, C.T)) |
Brian Silverman | e18cf50 | 2015-11-28 23:04:43 +0000 | [diff] [blame] | 142 | if controllability_rank != n: |
Austin Schuh | c9177b5 | 2015-11-28 13:18:31 -0800 | [diff] [blame] | 143 | glog.warning('Observability of %d != %d, unobservable state', |
Brian Silverman | e18cf50 | 2015-11-28 23:04:43 +0000 | [diff] [blame] | 144 | controllability_rank, n) |
Austin Schuh | e4a14f2 | 2015-03-01 00:12:29 -0800 | [diff] [blame] | 145 | |
Austin Schuh | 572ff40 | 2015-11-08 12:17:50 -0800 | [diff] [blame] | 146 | # Compute the steady state covariance matrix. |
Austin Schuh | c9177b5 | 2015-11-28 13:18:31 -0800 | [diff] [blame] | 147 | P_prior, rcond, w, S, T = slycot.sb02od(n=n, m=m, A=A.T, B=C.T, Q=Q, R=R, dico='D') |
Austin Schuh | 572ff40 | 2015-11-08 12:17:50 -0800 | [diff] [blame] | 148 | S = C * P_prior * C.T + R |
| 149 | K = numpy.linalg.lstsq(S.T, (P_prior * C.T).T)[0].T |
| 150 | P = (I - K * C) * P_prior |
Austin Schuh | e4a14f2 | 2015-03-01 00:12:29 -0800 | [diff] [blame] | 151 | |
| 152 | return K, P |
Austin Schuh | 86093ad | 2016-02-06 14:29:34 -0800 | [diff] [blame] | 153 | |
Austin Schuh | 3ad5ed8 | 2017-02-25 21:36:19 -0800 | [diff] [blame] | 154 | |
| 155 | def kalmd(A_continuous, B_continuous, Q_continuous, R_continuous, dt): |
| 156 | """Converts a continuous time kalman filter to discrete time. |
| 157 | |
| 158 | Args: |
| 159 | A_continuous: The A continuous matrix |
| 160 | B_continuous: the B continuous matrix |
| 161 | Q_continuous: The continuous cost matrix |
| 162 | R_continuous: The R continuous matrix |
| 163 | dt: Timestep |
| 164 | |
| 165 | The math for this is from: |
| 166 | https://www.mathworks.com/help/control/ref/kalmd.html |
| 167 | |
| 168 | Returns: |
| 169 | The discrete matrices of A, B, Q, and R. |
| 170 | """ |
| 171 | # TODO(austin): Verify that the dimensions make sense. |
| 172 | number_of_states = A_continuous.shape[0] |
| 173 | number_of_inputs = B_continuous.shape[1] |
| 174 | M = numpy.zeros((len(A_continuous) + number_of_inputs, |
| 175 | len(A_continuous) + number_of_inputs)) |
| 176 | M[0:number_of_states, 0:number_of_states] = A_continuous |
| 177 | M[0:number_of_states, number_of_states:] = B_continuous |
| 178 | M_exp = scipy.linalg.expm(M * dt) |
| 179 | A_discrete = M_exp[0:number_of_states, 0:number_of_states] |
| 180 | B_discrete = numpy.matrix(M_exp[0:number_of_states, number_of_states:]) |
| 181 | Q_continuous = (Q_continuous + Q_continuous.T) / 2.0 |
| 182 | R_continuous = (R_continuous + R_continuous.T) / 2.0 |
| 183 | M = numpy.concatenate((-A_continuous, Q_continuous), axis=1) |
| 184 | M = numpy.concatenate( |
| 185 | (M, numpy.concatenate((numpy.matrix( |
| 186 | numpy.zeros((number_of_states, number_of_states))), |
| 187 | numpy.transpose(A_continuous)), axis = 1)), axis = 0) |
| 188 | phi = numpy.matrix(scipy.linalg.expm(M*dt)) |
| 189 | phi12 = phi[0:number_of_states, number_of_states:(2*number_of_states)] |
| 190 | phi22 = phi[number_of_states:2*number_of_states, number_of_states:2*number_of_states] |
| 191 | Q_discrete = phi22.T * phi12 |
| 192 | Q_discrete = (Q_discrete + Q_discrete.T) / 2.0 |
| 193 | R_discrete = R_continuous / dt |
| 194 | return (A_discrete, B_discrete, Q_discrete, R_discrete) |
| 195 | |
| 196 | |
Austin Schuh | 86093ad | 2016-02-06 14:29:34 -0800 | [diff] [blame] | 197 | def TwoStateFeedForwards(B, Q): |
| 198 | """Computes the feed forwards constant for a 2 state controller. |
| 199 | |
| 200 | This will take the form U = Kff * (R(n + 1) - A * R(n)), where Kff is the |
| 201 | feed-forwards constant. It is important that Kff is *only* computed off |
| 202 | the goal and not the feed back terms. |
| 203 | |
| 204 | Args: |
| 205 | B: numpy.Matrix[num_states, num_inputs] The B matrix. |
| 206 | Q: numpy.Matrix[num_states, num_states] The Q (cost) matrix. |
| 207 | |
| 208 | Returns: |
| 209 | numpy.Matrix[num_inputs, num_states] |
| 210 | """ |
| 211 | |
| 212 | # We want to find the optimal U such that we minimize the tracking cost. |
| 213 | # This means that we want to minimize |
| 214 | # (B * U - (R(n+1) - A R(n)))^T * Q * (B * U - (R(n+1) - A R(n))) |
| 215 | # TODO(austin): This doesn't take into account the cost of U |
| 216 | |
| 217 | return numpy.linalg.inv(B.T * Q * B) * B.T * Q.T |