-
Notifications
You must be signed in to change notification settings - Fork 29
/
opLU.m
90 lines (79 loc) · 3.1 KB
/
opLU.m
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
classdef opLU < opFactorization
%opLU Operator representing the LU factorization of a
% matrix with optional iterative refinement. If the matrix
% is known to be symmetric, opLDL will be more efficient.
% If in addition, the matrix is known to be definite,
% opChol will be more efficient.
%
% opLU(A) creates an operator for multiplication by the (pseudo-)
% inverse of the matrix A implicitly represented by its LU
% factorization. Optionally, iterative refinement is performed.
% Note that A is an explicit matrix.
%
% The following attributes may be changed by the user:
% * nitref : the maximum number of iterative refinement steps (3)
% * itref_tol : iterative refinement tolerance (1.0e-8)
% * force_itref : force iterative refinement (false)
%
% See also lu.
%
% Dominique Orban <dominique.orban@gerad.ca>, 2014.
%
% Copyright 2009, Ewout van den Berg and Michael P. Friedlander
% See the file COPYING.txt for full copyright information.
% Use the command 'spot.gpl' to locate this file.
% http://www.cs.ubc.ca/labs/scl/spot
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Properties
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
properties( SetAccess = private )
P % Permutation operator
Q % Permutation operator
L % Lower triangular factor
U % Upper triangular factor
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Methods - Public
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
methods
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% opLU. Constructor
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function op = opLU(A)
if nargin ~= 1
error('Invalid number of arguments.');
end
[m,n] = size(A);
% Construct operator
op = op@opFactorization('LU', m, n);
B = A;
if ~issparse(A)
B = sparse(A);
end
op.A = opMatrix(B);
[op.L, op.U, p, q] = lu(B, 'vector');
op.P = opPermutation(p);
op.Q = opPermutation(q);
op.Ainv = op.Q * inv(op.U) * inv(op.L) * op.P;
op.cflag = ~isreal(A);
end % function opLU
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% transpose
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function opOut = transpose(op)
opOut = op.P' * inv(op.L.') * inv(op.U.') * op.Q';
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% conj
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function opOut = conj(op)
opOut = op.Q * inv(conj(op.U)) * inv(conj(op.L)) * op.P;
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ctranpose
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function opOut = ctranspose(op)
opOut = op.P' * inv(op.L') * inv(op.U') * op.Q';
end
end % methods - public
end % classdef