Mercurial > octave
annotate libinterp/corefcn/mgorth.cc @ 20802:8bb38ba1bad6
eliminate return statements after calls to print_usage
* __contourc__.cc, __dispatch__.cc, __dsearchn__.cc, __ichol__.cc,
__lin_interpn__.cc, __qp__.cc, balance.cc, betainc.cc, bsxfun.cc,
colloc.cc, daspk.cc, dasrt.cc, dassl.cc, defaults.cc, det.cc,
dlmread.cc, dot.cc, eig.cc, ellipj.cc, fft.cc, fft2.cc, fftn.cc,
filter.cc, find.cc, gcd.cc, givens.cc, hex2num.cc, inv.cc, lookup.cc,
lu.cc, max.cc, mgorth.cc, ordschur.cc, pinv.cc, profiler.cc, quad.cc,
qz.cc, rcond.cc, schur.cc, str2double.cc:
Eliminate return statements after calls to print_usage.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Fri, 04 Dec 2015 12:03:44 -0500 |
parents | e914b5399c67 |
children | f428cbe7576f |
rev | line source |
---|---|
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
1 /* |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
2 |
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
17787
diff
changeset
|
3 Copyright (C) 2009-2015 Carlo de Falco |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
4 Copyright (C) 2010 VZLU Prague |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
5 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
6 This file is part of Octave. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
7 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
8 Octave is free software; you can redistribute it and/or modify it |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
9 under the terms of the GNU General Public License as published by the |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
10 Free Software Foundation; either version 3 of the License, or (at your |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
11 option) any later version. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
12 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
13 Octave is distributed in the hope that it will be useful, but WITHOUT |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
16 for more details. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
17 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
18 You should have received a copy of the GNU General Public License |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
19 along with Octave; see the file COPYING. If not, see |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
20 <http://www.gnu.org/licenses/>. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
21 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
22 */ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
23 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
24 #ifdef HAVE_CONFIG_H |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
25 #include <config.h> |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
26 #endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
27 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
28 #include "oct-norm.h" |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14846
diff
changeset
|
29 #include "defun.h" |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
30 #include "error.h" |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
31 #include "gripes.h" |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
32 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
33 template <class ColumnVector, class Matrix, class RowVector> |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
12509
diff
changeset
|
34 static void |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
12509
diff
changeset
|
35 do_mgorth (ColumnVector& x, const Matrix& V, RowVector& h) |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
36 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
37 octave_idx_type Vc = V.columns (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
38 h = RowVector (Vc + 1); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
39 for (octave_idx_type j = 0; j < Vc; j++) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
40 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
41 ColumnVector Vcj = V.column (j); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
42 h(j) = RowVector (Vcj.hermitian ()) * x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
43 x -= h(j) * Vcj; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
44 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
45 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
46 h(Vc) = xnorm (x); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
47 if (real (h(Vc)) > 0) |
20230
e914b5399c67
Use in-place operators in C++ code where possible.
Rik <rik@octave.org>
parents:
20172
diff
changeset
|
48 x /= h(Vc); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
49 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
50 |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14846
diff
changeset
|
51 DEFUN (mgorth, args, nargout, |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
52 "-*- texinfo -*-\n\ |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14846
diff
changeset
|
53 @deftypefn {Built-in Function} {[@var{y}, @var{h}] =} mgorth (@var{x}, @var{v})\n\ |
15405
6d4cd994ea71
corrections to the doc for mgorth
Ed Meyer <eem2314@gmail.com>
parents:
15195
diff
changeset
|
54 Orthogonalize a given column vector @var{x} with respect to a set of\n\ |
20172
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
55 orthonormal vectors comprising the columns of @var{v} using the modified\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
56 Gram-Schmidt method.\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
57 \n\ |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
58 On exit, @var{y} is a unit vector such that:\n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
59 \n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
60 @example\n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
61 @group\n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
62 norm (@var{y}) = 1\n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
63 @var{v}' * @var{y} = 0\n\ |
15405
6d4cd994ea71
corrections to the doc for mgorth
Ed Meyer <eem2314@gmail.com>
parents:
15195
diff
changeset
|
64 @var{x} = [@var{v}, @var{y}]*@var{h}'\n\ |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
65 @end group\n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
66 @end example\n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
67 \n\ |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
68 @end deftypefn") |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
69 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
70 octave_value_list retval; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
71 |
14846
460a3c6d8bf1
maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents:
14501
diff
changeset
|
72 int nargin = args.length (); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
73 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
74 if (nargin != 2 || nargout > 2) |
20802
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20230
diff
changeset
|
75 print_usage (); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
76 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
77 octave_value arg_x = args(0); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
78 octave_value arg_v = args(1); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
79 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
80 if (arg_v.ndims () != 2 || arg_x.ndims () != 2 || arg_x.columns () != 1 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
81 || arg_v.rows () != arg_x.rows ()) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
82 { |
15405
6d4cd994ea71
corrections to the doc for mgorth
Ed Meyer <eem2314@gmail.com>
parents:
15195
diff
changeset
|
83 error ("mgorth: V should be a matrix, and X a column vector with" |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
84 " the same number of rows as V."); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
85 return retval; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
86 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
87 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
88 if (! arg_x.is_numeric_type () && ! arg_v.is_numeric_type ()) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
89 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
90 error ("mgorth: X and V must be numeric"); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
91 } |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
12509
diff
changeset
|
92 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
93 bool iscomplex = (arg_x.is_complex_type () || arg_v.is_complex_type ()); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
94 if (arg_x.is_single_type () || arg_v.is_single_type ()) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
95 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
96 if (iscomplex) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
97 { |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
98 FloatComplexColumnVector x |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
99 = arg_x.float_complex_column_vector_value (); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
100 FloatComplexMatrix V = arg_v.float_complex_matrix_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
101 FloatComplexRowVector h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
102 do_mgorth (x, V, h); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
103 retval(1) = h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
104 retval(0) = x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
105 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
106 else |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
107 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
108 FloatColumnVector x = arg_x.float_column_vector_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
109 FloatMatrix V = arg_v.float_matrix_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
110 FloatRowVector h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
111 do_mgorth (x, V, h); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
112 retval(1) = h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
113 retval(0) = x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
114 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
115 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
116 else |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
117 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
118 if (iscomplex) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
119 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
120 ComplexColumnVector x = arg_x.complex_column_vector_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
121 ComplexMatrix V = arg_v.complex_matrix_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
122 ComplexRowVector h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
123 do_mgorth (x, V, h); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
124 retval(1) = h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
125 retval(0) = x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
126 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
127 else |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
128 { |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
129 ColumnVector x = arg_x.column_vector_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
130 Matrix V = arg_v.matrix_value (); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
131 RowVector h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
132 do_mgorth (x, V, h); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
133 retval(1) = h; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
134 retval(0) = x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
135 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
136 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
137 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
138 return retval; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
139 } |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
140 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
141 /* |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
142 %!test |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
143 %! for ii=1:100 |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
144 %! assert (abs (mgorth (randn (5, 1), eye (5, 4))), [0 0 0 0 1]', eps); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
145 %! endfor |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
146 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
147 %!test |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
148 %! a = hilb (5); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
149 %! a(:, 1) /= norm (a(:, 1)); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
150 %! for ii = 1:5 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
151 %! a(:, ii) = mgorth (a(:, ii), a(:, 1:ii-1)); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
152 %! endfor |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
153 %! assert (a' * a, eye (5), 1e10); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
154 */ |