1 ## Copyright (C) 2009, 2010 Lukas F. Reichlin
3 ## This file is part of LTI Syncope.
5 ## LTI Syncope is free software: you can redistribute it and/or modify
6 ## it under the terms of the GNU General Public License as published by
7 ## the Free Software Foundation, either version 3 of the License, or
8 ## (at your option) any later version.
10 ## LTI Syncope is distributed in the hope that it will be useful,
11 ## but WITHOUT ANY WARRANTY; without even the implied warranty of
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 ## GNU General Public License for more details.
15 ## You should have received a copy of the GNU General Public License
16 ## along with LTI Syncope. If not, see <http://www.gnu.org/licenses/>.
19 ## @deftypefn {Function File} {[@var{bool}, @var{ncon}] =} isctrb (@var{sys})
20 ## @deftypefnx {Function File} {[@var{bool}, @var{ncon}] =} isctrb (@var{sys}, @var{tol})
21 ## @deftypefnx {Function File} {[@var{bool}, @var{ncon}] =} isctrb (@var{a}, @var{b})
22 ## @deftypefnx {Function File} {[@var{bool}, @var{ncon}] =} isctrb (@var{a}, @var{b}, @var{e})
23 ## @deftypefnx {Function File} {[@var{bool}, @var{ncon}] =} isctrb (@var{a}, @var{b}, @var{[]}, @var{tol})
24 ## @deftypefnx {Function File} {[@var{bool}, @var{ncon}] =} isctrb (@var{a}, @var{b}, @var{e}, @var{tol})
25 ## Logical check for system controllability.
26 ## For numerical reasons, @code{isctrb (sys)}
27 ## should be used instead of @code{rank (ctrb (sys))}.
32 ## LTI model. Descriptor state-space models are possible.
34 ## State transition matrix.
39 ## If @var{e} is empty @code{[]} or not specified, an identity matrix is assumed.
41 ## Optional roundoff parameter. Default value is 0.
47 ## System is not controllable.
49 ## System is controllable.
51 ## Number of controllable states.
54 ## @strong{Algorithm}@*
55 ## Uses SLICOT AB01OD and TG01HD by courtesy of
56 ## @uref{http://www.slicot.org, NICONET e.V.}
61 ## Author: Lukas Reichlin <lukas.reichlin@gmail.com>
62 ## Created: October 2009
65 function [bool, ncont] = isctrb (a, b = [], e = [], tol = [])
67 if (nargin < 1 || nargin > 4)
69 elseif (isa (a, "lti")) # isctrb (sys), isctrb (sys, tol)
74 [a, b, c, d, e] = dssdata (a, []);
75 elseif (nargin < 2) # isctrb (a, b), isctrb (a, b, tol)
77 elseif (! is_real_square_matrix (a) || ! is_real_matrix (b) || rows (a) != rows (b))
78 error ("isctrb: a(%dx%d), b(%dx%d) not conformal",
79 rows (a), columns (a), rows (b), columns (b));
80 elseif (! isempty (e) && (! is_real_square_matrix (e) || ! size_equal (e, a)))
81 error ("isctrb: a(%dx%d), e(%dx%d) not conformal",
82 rows (a), columns (a), rows (e), columns (e));
86 tol = 0; # default tolerance
87 elseif (! is_real_scalar (tol))
88 error ("isctrb: tol must be a real scalar");
92 [~, ~, ~, ncont] = slab01od (a, b, tol);
94 [~, ~, ~, ~, ~, ~, ncont] = sltg01hd (a, e, b, zeros (1, columns (a)), tol);
97 bool = (ncont == rows (a));