GNU Octave 7.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
__eigs__.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 2005-2022 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26#if defined (HAVE_CONFIG_H)
27# include "config.h"
28#endif
29
30#include <limits>
31#include <string>
32
33#include "Matrix.h"
34#include "eigs-base.h"
35#include "unwind-prot.h"
36
37#include "defun.h"
38#include "error.h"
39#include "errwarn.h"
40#include "interpreter-private.h"
41#include "oct-map.h"
42#include "ov.h"
43#include "ovl.h"
44#include "pager.h"
45#include "parse.h"
46#include "variables.h"
47
48OCTAVE_NAMESPACE_BEGIN
49
50#if defined (HAVE_ARPACK)
51
53{
54public:
55
57 eigs_func (const ColumnVector& x, int& eigs_error);
58
60 eigs_complex_func (const ComplexColumnVector& x, int& eigs_error);
61
62 //--------
63
64 // Pointer for user defined function.
66
67 // Have we warned about imaginary values returned from user function?
68 bool m_warned_imaginary = false;
69};
70
71// Is this a recursive call?
72static int call_depth = 0;
73
75eigs_callback::eigs_func (const ColumnVector& x, int& eigs_error)
76{
77 ColumnVector retval;
79 args(0) = x;
80
82 {
84
85 try
86 {
87 tmp = octave::feval (m_eigs_fcn, args, 1);
88 }
89 catch (octave::execution_exception& ee)
90 {
91 err_user_supplied_eval (ee, "eigs");
92 }
93
94 if (tmp.length () && tmp(0).is_defined ())
95 {
96 if (! m_warned_imaginary && tmp(0).iscomplex ())
97 {
98 warning ("eigs: ignoring imaginary part returned from user-supplied function");
99 m_warned_imaginary = true;
100 }
101
102 retval = tmp(0).xvector_value ("eigs: evaluation of user-supplied function failed");
103 }
104 else
105 {
106 eigs_error = 1;
107 err_user_supplied_eval ("eigs");
108 }
109 }
110
111 return retval;
112}
113
116 int& eigs_error)
117{
118 ComplexColumnVector retval;
120 args(0) = x;
121
122 if (m_eigs_fcn.is_defined ())
123 {
125
126 try
127 {
128 tmp = octave::feval (m_eigs_fcn, args, 1);
129 }
130 catch (octave::execution_exception& ee)
131 {
132 err_user_supplied_eval (ee, "eigs");
133 }
134
135 if (tmp.length () && tmp(0).is_defined ())
136 {
137 retval = tmp(0).xcomplex_vector_value ("eigs: evaluation of user-supplied function failed");
138 }
139 else
140 {
141 eigs_error = 1;
142 err_user_supplied_eval ("eigs");
143 }
144 }
145
146 return retval;
147}
148
149#endif
150
151DEFMETHOD (__eigs__, interp, args, nargout,
152 doc: /* -*- texinfo -*-
153@deftypefn {} {@var{d} =} __eigs__ (@var{A})
154@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{k})
155@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{k}, @var{sigma})
156@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{k}, @var{sigma}, @var{opts})
157@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B})
158@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B}, @var{k})
159@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B}, @var{k}, @var{sigma})
160@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B}, @var{k}, @var{sigma}, @var{opts})
161@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n})
162@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B})
163@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{k})
164@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B}, @var{k})
165@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{k}, @var{sigma})
166@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B}, @var{k}, @var{sigma})
167@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{k}, @var{sigma}, @var{opts})
168@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B}, @var{k}, @var{sigma}, @var{opts})
169@deftypefnx {} {[@var{V}, @var{d}] =} __eigs__ (@var{A}, @dots{})
170@deftypefnx {} {[@var{V}, @var{d}] =} __eigs__ (@var{af}, @var{n}, @dots{})
171@deftypefnx {} {[@var{V}, @var{d}, @var{flag}] =} __eigs__ (@var{A}, @dots{})
172@deftypefnx {} {[@var{V}, @var{d}, @var{flag}] =} __eigs__ (@var{af}, @var{n}, @dots{})
173Undocumented internal function.
174@end deftypefn */)
175{
176#if defined (HAVE_ARPACK)
177
178 int nargin = args.length ();
179
180 if (nargin == 0)
181 print_usage ();
182
183 octave_value_list retval;
184
185 std::string fcn_name;
186 octave_idx_type n = 0;
187 octave_idx_type k = 6;
188 Complex sigma = 0.0;
189 double sigmar, sigmai;
190 bool have_sigma = false;
191 std::string typ = "LM";
192 Matrix amm, bmm, bmt;
193 ComplexMatrix acm, bcm, bct;
194 SparseMatrix asmm, bsmm, bsmt;
195 SparseComplexMatrix ascm, bscm, bsct;
196 int b_arg = 0;
197 bool have_b = false;
198 bool have_a_fun = false;
199 bool a_is_complex = false;
200 bool b_is_complex = false;
201 bool symmetric = false;
202 bool sym_tested = false;
203 bool cholB = false;
204 bool a_is_sparse = false;
205 bool b_is_sparse = false;
206 ColumnVector permB;
207 int arg_offset = 0;
208 double tol = std::numeric_limits<double>::epsilon ();
209 int maxit = 300;
210 int disp = 0;
211 octave_idx_type p = -1;
212 ColumnVector resid;
213 ComplexColumnVector cresid;
214 octave_idx_type info = 1;
215
216 eigs_callback callback;
217
218 unwind_protect_var<int> restore_var (call_depth);
219 call_depth++;
220
221 if (call_depth > 1)
222 error ("eigs: invalid recursive call");
223
224 if (args(0).is_function_handle () || args(0).is_inline_function ()
225 || args(0).is_string ())
226 {
227 callback.m_eigs_fcn = get_function_handle (interp, args(0), "x");
228
229 if (callback.m_eigs_fcn.is_undefined ())
230 error ("eigs: unknown function");
231
232 if (nargin < 2)
233 error ("eigs: incorrect number of arguments");
234
235 n = args(1).nint_value ();
236 arg_offset = 1;
237 have_a_fun = true;
238 }
239 else
240 {
241 if (args(0).iscomplex ())
242 {
243 if (args(0).issparse ())
244 {
245 ascm = (args(0).sparse_complex_matrix_value ());
246 a_is_sparse = true;
247 }
248 else
249 acm = (args(0).complex_matrix_value ());
250 a_is_complex = true;
251 }
252 else
253 {
254 if (args(0).issparse ())
255 {
256 asmm = (args(0).sparse_matrix_value ());
257 a_is_sparse = true;
258 }
259 else
260 {
261 amm = (args(0).matrix_value ());
262 }
263 }
264 }
265
266 // Note hold off reading B until later to avoid issues of double
267 // copies of the matrix if B is full/real while A is complex.
268 if (nargin > 1 + arg_offset
269 && ! (args(1 + arg_offset).is_real_scalar ()))
270 {
271 if (args(1+arg_offset).iscomplex ())
272 {
273 b_arg = 1+arg_offset;
274 if (args(b_arg).issparse ())
275 {
276 bscm = (args(b_arg).sparse_complex_matrix_value ());
277 b_is_sparse = true;
278 }
279 else
280 bcm = (args(b_arg).complex_matrix_value ());
281 have_b = true;
282 b_is_complex = true;
283 arg_offset++;
284 }
285 else
286 {
287 b_arg = 1+arg_offset;
288 if (args(b_arg).issparse ())
289 {
290 bsmm = (args(b_arg).sparse_matrix_value ());
291 b_is_sparse = true;
292 }
293 else
294 bmm = (args(b_arg).matrix_value ());
295 have_b = true;
296 arg_offset++;
297 }
298 }
299
300 if (nargin > (1+arg_offset))
301 k = args(1+arg_offset).nint_value ();
302
303 if (nargin > (2+arg_offset))
304 {
305 if (args(2+arg_offset).is_string ())
306 {
307 typ = args(2+arg_offset).string_value ();
308
309 // Use STL function to convert to upper case
310 transform (typ.begin (), typ.end (), typ.begin (), toupper);
311
312 sigma = 0.0;
313 }
314 else
315 {
316 sigma = args(2+arg_offset).xcomplex_value ("eigs: SIGMA must be a scalar or a string");
317
318 have_sigma = true;
319 }
320 }
321
322 sigmar = sigma.real ();
323 sigmai = sigma.imag ();
324
325 if (nargin > (3+arg_offset))
326 {
327 if (! args(3+arg_offset).isstruct ())
328 error ("eigs: OPTS argument must be a structure");
329
330 octave_scalar_map map = args(3+arg_offset).xscalar_map_value ("eigs: OPTS argument must be a scalar structure");
331
332 octave_value tmp;
333
334 // issym is ignored for complex matrix inputs
335 tmp = map.getfield ("issym");
336 if (tmp.is_defined ())
337 {
338 if (tmp.numel () != 1)
339 error ("eigs: OPTS.issym must be a scalar value");
340
341 symmetric = tmp.xbool_value ("eigs: OPTS.issym must be a logical value");
342 sym_tested = true;
343 }
344
345 // isreal is ignored if A is not a function
346 if (have_a_fun)
347 {
348 tmp = map.getfield ("isreal");
349 if (tmp.is_defined ())
350 {
351 if (tmp.numel () != 1)
352 error ("eigs: OPTS.isreal must be a scalar value");
353
354 a_is_complex = ! tmp.xbool_value ("eigs: OPTS.isreal must be a logical value");
355 }
356 }
357
358 tmp = map.getfield ("tol");
359 if (tmp.is_defined ())
360 tol = tmp.double_value ();
361
362 tmp = map.getfield ("maxit");
363 if (tmp.is_defined ())
364 maxit = tmp.nint_value ();
365
366 tmp = map.getfield ("p");
367 if (tmp.is_defined ())
368 p = tmp.nint_value ();
369
370 tmp = map.getfield ("v0");
371 if (tmp.is_defined ())
372 {
373 if (a_is_complex || b_is_complex)
375 else
376 resid = ColumnVector (tmp.vector_value ());
377 }
378
379 tmp = map.getfield ("disp");
380 if (tmp.is_defined ())
381 disp = tmp.nint_value ();
382
383 tmp = map.getfield ("cholB");
384 if (tmp.is_defined ())
385 {
386 if (tmp.numel () != 1)
387 error ("eigs: OPTS.cholB must be a scalar value");
388
389 cholB = tmp.xbool_value ("eigs: OPTS.cholB must be a logical value");
390 }
391
392 tmp = map.getfield ("permB");
393 if (tmp.is_defined ())
394 permB = ColumnVector (tmp.vector_value ()) - 1.0;
395 }
396
397 if (nargin > (4+arg_offset))
398 error ("eigs: incorrect number of arguments");
399
400 // Test undeclared (no issym) matrix inputs for symmetry
401 if (! sym_tested && ! have_a_fun)
402 {
403 if (a_is_complex)
404 {
405 if (a_is_sparse)
406 symmetric = ascm.ishermitian ();
407 else
408 symmetric = acm.ishermitian ();
409 }
410 else
411 {
412 if (a_is_sparse)
413 symmetric = asmm.issymmetric ();
414 else
415 symmetric = amm.issymmetric ();
416 }
417 }
418
419 if (have_b)
420 {
421 if (a_is_complex || b_is_complex)
422 {
423 if (b_is_sparse)
424 bscm = args(b_arg).sparse_complex_matrix_value ();
425 else
426 bcm = args(b_arg).complex_matrix_value ();
427 }
428 else
429 {
430 if (b_is_sparse)
431 bsmm = args(b_arg).sparse_matrix_value ();
432 else
433 bmm = args(b_arg).matrix_value ();
434 }
435 }
436
437 // Mode 1 for SM mode seems unstable for some reason.
438 // Use Mode 3 instead, with sigma = 0.
439 if (! have_sigma && typ == "SM")
440 have_sigma = true;
441
442 octave_idx_type nconv;
443 if (a_is_complex || b_is_complex)
444 {
446 eigs_complex_fcn = [&callback] (const ComplexColumnVector& x,
447 int& eigs_error)
448 {
449 return callback.eigs_complex_func (x, eigs_error);
450 };
451
452 ComplexMatrix eig_vec;
453 ComplexColumnVector eig_val;
454
455 if (have_a_fun)
456 {
457 if (b_is_sparse)
459 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
460 eig_val, bscm, permB, cresid, octave_stdout, tol,
461 (nargout > 1), cholB, disp, maxit);
462 else
464 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
465 eig_val, bcm, permB, cresid, octave_stdout, tol,
466 (nargout > 1), cholB, disp, maxit);
467 }
468 else if (have_sigma)
469 {
470 if (a_is_sparse)
472 (ascm, sigma, k, p, info, eig_vec, eig_val, bscm, permB,
473 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
474 maxit);
475 else
477 (acm, sigma, k, p, info, eig_vec, eig_val, bcm, permB,
478 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
479 maxit);
480 }
481 else
482 {
483 if (a_is_sparse)
485 (ascm, typ, k, p, info, eig_vec, eig_val, bscm, permB,
486 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
487 maxit);
488 else
490 (acm, typ, k, p, info, eig_vec, eig_val, bcm, permB,
491 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
492 maxit);
493 }
494
495 if (nargout < 2)
496 {
497 if (symmetric)
498 retval(0) = real (eig_val);
499 else
500 retval(0) = eig_val;
501 }
502 else
503 {
504 if (symmetric)
505 retval = ovl (eig_vec, DiagMatrix (real (eig_val)), double (info));
506 else
507 retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
508 }
509 }
510 else if (sigmai != 0.0)
511 {
513 eigs_complex_fcn = [&callback] (const ComplexColumnVector& x,
514 int& eigs_error)
515 {
516 return callback.eigs_complex_func (x, eigs_error);
517 };
518
519 // Promote real problem to a complex one.
520 ComplexMatrix eig_vec;
521 ComplexColumnVector eig_val;
522
523 if (have_a_fun)
524 {
525 if (b_is_sparse)
527 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
528 eig_val, bscm, permB, cresid, octave_stdout, tol,
529 (nargout > 1), cholB, disp, maxit);
530 else
532 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
533 eig_val, bcm, permB, cresid, octave_stdout, tol,
534 (nargout > 1), cholB, disp, maxit);
535 }
536 else
537 {
538 if (a_is_sparse)
540 (SparseComplexMatrix (asmm), sigma, k, p, info, eig_vec,
541 eig_val, SparseComplexMatrix (bsmm), permB, cresid,
542 octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
543 else
545 (ComplexMatrix (amm), sigma, k, p, info, eig_vec,
546 eig_val, ComplexMatrix (bmm), permB, cresid,
547 octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
548 }
549
550 if (nargout < 2)
551 {
552 if (symmetric)
553 retval(0) = real (eig_val);
554 else
555 retval(0) = eig_val;
556 }
557 else
558 {
559 if (symmetric)
560 retval = ovl (eig_vec, DiagMatrix (real (eig_val)), double (info));
561 else
562 retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
563 }
564 }
565 else
566 {
567 EigsFunc eigs_fcn = [&callback] (const ColumnVector& x, int& eigs_error)
568 {
569 return callback.eigs_func (x, eigs_error);
570 };
571
572 if (symmetric)
573 {
574 Matrix eig_vec;
575 ColumnVector eig_val;
576
577 if (have_a_fun)
578 {
579 if (b_is_sparse)
581 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
582 eig_val, bsmm, permB, resid, octave_stdout, tol,
583 (nargout > 1), cholB, disp, maxit);
584 else
586 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
587 eig_val, bmm, permB, resid, octave_stdout, tol,
588 (nargout > 1), cholB, disp, maxit);
589 }
590 else if (have_sigma)
591 {
592 if (a_is_sparse)
594 (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm,
595 permB, resid, octave_stdout, tol, (nargout > 1),
596 cholB, disp, maxit);
597 else
599 (amm, sigmar, k, p, info, eig_vec, eig_val, bmm,
600 permB, resid, octave_stdout, tol, (nargout > 1),
601 cholB, disp, maxit);
602 }
603 else
604 {
605 if (a_is_sparse)
607 (asmm, typ, k, p, info, eig_vec, eig_val, bsmm,
608 permB, resid, octave_stdout, tol, (nargout > 1),
609 cholB, disp, maxit);
610 else
612 (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
613 resid, octave_stdout, tol, (nargout > 1), cholB,
614 disp, maxit);
615 }
616
617 if (nargout < 2)
618 retval(0) = eig_val;
619 else
620 retval = ovl (eig_vec, DiagMatrix (eig_val), double (info));
621 }
622 else
623 {
624 ComplexMatrix eig_vec;
625 ComplexColumnVector eig_val;
626
627 if (have_a_fun)
628 {
629 if (b_is_sparse)
631 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
632 eig_val, bsmm, permB, resid, octave_stdout, tol,
633 (nargout > 1), cholB, disp, maxit);
634 else
636 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
637 eig_val, bmm, permB, resid, octave_stdout, tol,
638 (nargout > 1), cholB, disp, maxit);
639 }
640 else if (have_sigma)
641 {
642 if (a_is_sparse)
644 (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm,
645 permB, resid, octave_stdout, tol, (nargout > 1),
646 cholB, disp, maxit);
647 else
649 (amm, sigmar, k, p, info, eig_vec, eig_val, bmm,
650 permB, resid, octave_stdout, tol, (nargout > 1),
651 cholB, disp, maxit);
652 }
653 else
654 {
655 if (a_is_sparse)
657 (asmm, typ, k, p, info, eig_vec, eig_val, bsmm,
658 permB, resid, octave_stdout, tol, (nargout > 1),
659 cholB, disp, maxit);
660 else
662 (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
663 resid, octave_stdout, tol, (nargout > 1), cholB,
664 disp, maxit);
665 }
666
667 if (nargout < 2)
668 retval(0) = eig_val;
669 else
670 retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
671 }
672 }
673
674 if (nconv <= 0)
675 warning_with_id ("Octave:eigs:UnconvergedEigenvalues",
676 "eigs: None of the %" OCTAVE_IDX_TYPE_FORMAT
677 " requested eigenvalues converged", k);
678 else if (nconv < k)
679 warning_with_id ("Octave:eigs:UnconvergedEigenvalues",
680 "eigs: Only %" OCTAVE_IDX_TYPE_FORMAT
681 " of the %" OCTAVE_IDX_TYPE_FORMAT
682 " requested eigenvalues converged",
683 nconv, k);
684
685 if (! fcn_name.empty ())
686 {
687 symbol_table& symtab = interp.get_symbol_table ();
688
689 symtab.clear_function (fcn_name);
690 }
691
692 return retval;
693
694#else
695
696 octave_unused_parameter (interp);
697 octave_unused_parameter (args);
698 octave_unused_parameter (nargout);
699
700 err_disabled_feature ("eigs", "ARPACK");
701
702#endif
703}
704
705/*
706## No test needed for internal helper function.
707%!assert (1)
708*/
709
710OCTAVE_NAMESPACE_END
static int call_depth
Definition: __eigs__.cc:72
OCTAVE_API bool ishermitian(void) const
Definition: CMatrix.cc:174
Definition: dMatrix.h:42
OCTAVE_API bool issymmetric(void) const
Definition: dMatrix.cc:136
OCTAVE_API bool ishermitian(void) const
Definition: CSparse.cc:142
OCTAVE_API bool issymmetric(void) const
Definition: dSparse.cc:131
octave_value getfield(const std::string &key) const
Definition: oct-map.cc:183
octave_idx_type length(void) const
Definition: ovl.h:113
OCTINTERP_API Array< double > vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
octave_idx_type numel(void) const
Definition: ov.h:604
bool is_defined(void) const
Definition: ov.h:637
int nint_value(bool frc_str_conv=false) const
Definition: ov.h:864
OCTINTERP_API Array< Complex > complex_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
OCTINTERP_API bool xbool_value(const char *fmt,...) const
bool is_undefined(void) const
Definition: ov.h:640
double double_value(bool frc_str_conv=false) const
Definition: ov.h:886
void clear_function(const std::string &name)
Definition: symtab.cc:434
ColumnVector real(const ComplexColumnVector &a)
Definition: dColVector.cc:137
OCTINTERP_API void print_usage(void)
Definition: defun-int.h:72
#define DEFMETHOD(name, interp_name, args_name, nargout_name, doc)
Macro to define a builtin method.
Definition: defun.h:111
octave_idx_type EigsRealSymmetricFunc(EigsFunc fun, octave_idx_type n_arg, const std::string &_typ, double sigma, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, Matrix &eig_vec, ColumnVector &eig_val, const M &_b, ColumnVector &permB, ColumnVector &resid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:1381
octave_idx_type EigsRealNonSymmetricMatrixShift(const M &m, double sigmar, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, ComplexMatrix &eig_vec, ComplexColumnVector &eig_val, const M &_b, ColumnVector &permB, ColumnVector &resid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:2113
octave_idx_type EigsRealSymmetricMatrixShift(const M &m, double sigma, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, Matrix &eig_vec, ColumnVector &eig_val, const M &_b, ColumnVector &permB, ColumnVector &resid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:1078
octave_idx_type EigsRealNonSymmetricMatrix(const M &m, const std::string typ, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, ComplexMatrix &eig_vec, ComplexColumnVector &eig_val, const M &_b, ColumnVector &permB, ColumnVector &resid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:1762
octave_idx_type EigsComplexNonSymmetricMatrix(const M &m, const std::string typ, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, ComplexMatrix &eig_vec, ComplexColumnVector &eig_val, const M &_b, ColumnVector &permB, ComplexColumnVector &cresid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:2927
octave_idx_type EigsComplexNonSymmetricFunc(EigsComplexFunc fun, octave_idx_type n_arg, const std::string &_typ, Complex sigma, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, ComplexMatrix &eig_vec, ComplexColumnVector &eig_val, const M &_b, ColumnVector &permB, ComplexColumnVector &cresid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:3553
octave_idx_type EigsRealNonSymmetricFunc(EigsFunc fun, octave_idx_type n_arg, const std::string &_typ, double sigmar, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, ComplexMatrix &eig_vec, ComplexColumnVector &eig_val, const M &_b, ColumnVector &permB, ColumnVector &resid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:2485
octave_idx_type EigsComplexNonSymmetricMatrixShift(const M &m, Complex sigma, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, ComplexMatrix &eig_vec, ComplexColumnVector &eig_val, const M &_b, ColumnVector &permB, ComplexColumnVector &cresid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:3230
octave_idx_type EigsRealSymmetricMatrix(const M &m, const std::string typ, octave_idx_type k_arg, octave_idx_type p_arg, octave_idx_type &info, Matrix &eig_vec, ColumnVector &eig_val, const M &_b, ColumnVector &permB, ColumnVector &resid, std::ostream &os, double tol, bool rvec, bool cholB, int disp, int maxit)
Definition: eigs-base.cc:786
std::function< ColumnVector(const ColumnVector &x, int &eigs_error)> EigsFunc
Definition: eigs-base.h:39
std::function< ComplexColumnVector(const ComplexColumnVector &x, int &eigs_error)> EigsComplexFunc
Definition: eigs-base.h:43
void warning(const char *fmt,...)
Definition: error.cc:1055
void warning_with_id(const char *id, const char *fmt,...)
Definition: error.cc:1070
void error(const char *fmt,...)
Definition: error.cc:980
void err_disabled_feature(const std::string &fcn, const std::string &feature, const std::string &pkg)
Definition: errwarn.cc:53
void err_user_supplied_eval(const char *name)
Definition: errwarn.cc:152
ColumnVector transform(const Matrix &m, double x, double y, double z)
Definition: graphics.cc:5861
F77_RET_T const F77_DBLE * x
class OCTAVE_API ComplexColumnVector
Definition: mx-fwd.h:46
class OCTAVE_API ColumnVector
Definition: mx-fwd.h:45
OCTINTERP_API octave_value_list feval(const char *name, const octave_value_list &args=octave_value_list(), int nargout=0)
octave_value get_function_handle(interpreter &interp, const octave_value &arg, const std::string &parameter_name)
std::complex< double > Complex
Definition: oct-cmplx.h:33
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition: ovl.h:211
#define octave_stdout
Definition: pager.h:314
bool m_warned_imaginary
Definition: __eigs__.cc:68
ComplexColumnVector eigs_complex_func(const ComplexColumnVector &x, int &eigs_error)
Definition: __eigs__.cc:115
ColumnVector eigs_func(const ColumnVector &x, int &eigs_error)
Definition: __eigs__.cc:75
octave_value m_eigs_fcn
Definition: __eigs__.cc:65