GNU Octave 11.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
__eigs__.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 2005-2026 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
49
50#if defined (HAVE_ARPACK)
51
52struct eigs_callback
53{
54public:
55
56 eigs_callback (octave::interpreter& interp)
57 : m_interpreter (interp)
58 { }
59
61 eigs_func (const ColumnVector& x, int& eigs_error);
62
64 eigs_complex_func (const ComplexColumnVector& x, int& eigs_error);
65
66 //--------
67
68 octave::interpreter& m_interpreter;
69
70 // Pointer for user defined function.
71 octave_value m_eigs_fcn;
72
73 // Have we warned about imaginary values returned from user function?
74 bool m_warned_imaginary = false;
75};
76
77// Is this a recursive call?
78static int call_depth = 0;
79
81eigs_callback::eigs_func (const ColumnVector& x, int& eigs_error)
82{
83 ColumnVector retval;
85 args(0) = x;
86
87 if (m_eigs_fcn.is_defined ())
88 {
90
91 try
92 {
93 tmp = m_interpreter.feval (m_eigs_fcn, args, 1);
94 }
95 catch (octave::execution_exception& ee)
96 {
97 err_user_supplied_eval (ee, "eigs");
98 }
99
100 if (tmp.length () && tmp(0).is_defined ())
101 {
102 if (! m_warned_imaginary && tmp(0).iscomplex ())
103 {
104 warning ("eigs: ignoring imaginary part returned from user-supplied function");
105 m_warned_imaginary = true;
106 }
107
108 retval = tmp(0).xvector_value ("eigs: evaluation of user-supplied function failed");
109 }
110 else
111 {
112 eigs_error = 1;
113 err_user_supplied_eval ("eigs");
114 }
115 }
116
117 return retval;
118}
119
121eigs_callback::eigs_complex_func (const ComplexColumnVector& x,
122 int& eigs_error)
123{
124 ComplexColumnVector retval;
126 args(0) = x;
127
128 if (m_eigs_fcn.is_defined ())
129 {
131
132 try
133 {
134 tmp = m_interpreter.feval (m_eigs_fcn, args, 1);
135 }
136 catch (octave::execution_exception& ee)
137 {
138 err_user_supplied_eval (ee, "eigs");
139 }
140
141 if (tmp.length () && tmp(0).is_defined ())
142 {
143 retval = tmp(0).xcomplex_vector_value ("eigs: evaluation of user-supplied function failed");
144 }
145 else
146 {
147 eigs_error = 1;
148 err_user_supplied_eval ("eigs");
149 }
150 }
151
152 return retval;
153}
154
155#endif
156
157DEFMETHOD (__eigs__, interp, args, nargout,
158 doc: /* -*- texinfo -*-
159@deftypefn {} {@var{d} =} __eigs__ (@var{A})
160@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{k})
161@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{k}, @var{sigma})
162@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{k}, @var{sigma}, @var{opts})
163@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B})
164@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B}, @var{k})
165@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B}, @var{k}, @var{sigma})
166@deftypefnx {} {@var{d} =} __eigs__ (@var{A}, @var{B}, @var{k}, @var{sigma}, @var{opts})
167@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n})
168@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B})
169@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{k})
170@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B}, @var{k})
171@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{k}, @var{sigma})
172@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B}, @var{k}, @var{sigma})
173@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{k}, @var{sigma}, @var{opts})
174@deftypefnx {} {@var{d} =} __eigs__ (@var{af}, @var{n}, @var{B}, @var{k}, @var{sigma}, @var{opts})
175@deftypefnx {} {[@var{V}, @var{d}] =} __eigs__ (@var{A}, @dots{})
176@deftypefnx {} {[@var{V}, @var{d}] =} __eigs__ (@var{af}, @var{n}, @dots{})
177@deftypefnx {} {[@var{V}, @var{d}, @var{flag}] =} __eigs__ (@var{A}, @dots{})
178@deftypefnx {} {[@var{V}, @var{d}, @var{flag}] =} __eigs__ (@var{af}, @var{n}, @dots{})
179Undocumented internal function.
180@end deftypefn */)
181{
182#if defined (HAVE_ARPACK)
183
184 int nargin = args.length ();
185
186 if (nargin == 0)
187 print_usage ();
188
189 octave_value_list retval;
190
191 std::string fcn_name;
192 octave_idx_type n = 0;
193 octave_idx_type k = 6;
194 Complex sigma = 0.0;
195 double sigmar, sigmai;
196 bool have_sigma = false;
197 std::string typ = "LM";
198 Matrix amm, bmm, bmt;
199 ComplexMatrix acm, bcm, bct;
200 SparseMatrix asmm, bsmm, bsmt;
201 SparseComplexMatrix ascm, bscm, bsct;
202 int b_arg = 0;
203 bool have_b = false;
204 bool have_a_fcn = false;
205 bool a_is_complex = false;
206 bool b_is_complex = false;
207 bool symmetric = false;
208 bool sym_tested = false;
209 bool cholB = false;
210 bool a_is_sparse = false;
211 bool b_is_sparse = false;
212 ColumnVector permB;
213 int arg_offset = 0;
214 double tol = std::numeric_limits<double>::epsilon ();
215 int maxit = 300;
216 int disp = 0;
217 octave_idx_type p = -1;
218 ColumnVector resid;
219 ComplexColumnVector cresid;
220 octave_idx_type info = 1;
221
222 eigs_callback callback (interp);
223
224 unwind_protect_var<int> restore_var (call_depth);
225 call_depth++;
226
227 if (call_depth > 1)
228 error ("eigs: invalid recursive call");
229
230 if (args(0).is_function_handle () || args(0).is_inline_function ()
231 || args(0).is_string ())
232 {
233 have_a_fcn = true;
234 callback.m_eigs_fcn = get_function_handle (interp, args(0), "x");
235
236 if (callback.m_eigs_fcn.is_undefined ())
237 error ("eigs: unknown function");
238 if (nargin < 2)
239 error ("eigs: incorrect number of arguments");
240
241 n = args(1).strict_int_value ();
242 arg_offset = 1;
243 }
244 else
245 {
246 // Matrix A
247 if (args(0).iscomplex ())
248 {
249 a_is_complex = true;
250 if (args(0).issparse ())
251 {
252 a_is_sparse = true;
253 ascm = args(0).sparse_complex_matrix_value ();
254 // Validation is 17x faster in C++ than in m-file.
255 if (ascm.any_element_is_inf_or_nan ())
256 error ("eigs: matrix A contains Inf or NaN values");
257 }
258 else
259 {
260 acm = args(0).complex_matrix_value ();
261 if (acm.any_element_is_inf_or_nan ())
262 error ("eigs: matrix A contains Inf or NaN values");
263 }
264 }
265 else
266 {
267 if (args(0).issparse ())
268 {
269 a_is_sparse = true;
270 asmm = args(0).sparse_matrix_value ();
271 if (asmm.any_element_is_inf_or_nan ())
272 error ("eigs: matrix A contains Inf or NaN values");
273 }
274 else
275 {
276 amm = args(0).matrix_value ();
277 if (amm.any_element_is_inf_or_nan ())
278 error ("eigs: matrix A contains Inf or NaN values");
279 }
280 }
281 }
282
283 // Note: hold off converting B until later to avoid issues of double
284 // copies of the matrix if B is full/real while A is complex.
285 if (nargin > (1+arg_offset) && ! args(1+arg_offset).is_real_scalar ())
286 {
287 have_b = true;
288 b_arg = 1 + arg_offset;
289 if (args(b_arg).iscomplex ())
290 {
291 b_is_complex = true;
292 if (args(b_arg).issparse ())
293 b_is_sparse = true;
294 }
295 else
296 {
297 if (args(b_arg).issparse ())
298 b_is_sparse = true;
299 }
300 arg_offset++;
301 }
302
303 if (nargin > (1+arg_offset))
304 k = args(1+arg_offset).strict_int_value ();
305
306 if (nargin > (2+arg_offset))
307 {
308 if (args(2+arg_offset).is_string ())
309 {
310 typ = args(2+arg_offset).string_value ();
311 // Use STL function to convert to upper case
312 transform (typ.begin (), typ.end (), typ.begin (), toupper);
313
314 sigma = 0.0;
315 }
316 else
317 {
318 have_sigma = true;
319 sigma = args(2+arg_offset).xcomplex_value ("eigs: SIGMA must be a scalar or a string");
320 }
321 }
322
323 sigmar = sigma.real ();
324 sigmai = sigma.imag ();
325
326 if (nargin > (3+arg_offset))
327 {
328 if (! args(3+arg_offset).isstruct ())
329 error ("eigs: OPTS argument must be a structure");
330
331 octave_scalar_map map = args(3+arg_offset).xscalar_map_value ("eigs: OPTS argument must be a scalar structure");
332
333 octave_value tmp;
334
335 // issym is ignored for complex matrix inputs
336 tmp = map.getfield ("issym");
337 if (tmp.is_defined ())
338 {
339 if (tmp.numel () != 1)
340 error ("eigs: OPTS.issym must be a scalar value");
341
342 symmetric = tmp.strict_bool_value ("eigs: OPTS.issym must be a logical value");
343 sym_tested = true;
344 }
345
346 // isreal is ignored if A is not a function
347 if (have_a_fcn)
348 {
349 tmp = map.getfield ("isreal");
350 if (tmp.is_defined ())
351 {
352 if (tmp.numel () != 1)
353 error ("eigs: OPTS.isreal must be a scalar value");
354
355 a_is_complex = ! tmp.strict_bool_value ("eigs: OPTS.isreal must be a logical value");
356 }
357 }
358
359 tmp = map.getfield ("tol");
360 if (tmp.is_defined ())
361 tol = tmp.double_value ();
362
363 tmp = map.getfield ("maxit");
364 if (tmp.is_defined ())
365 maxit = tmp.strict_int_value ();
366
367 tmp = map.getfield ("p");
368 if (tmp.is_defined ())
369 p = tmp.strict_int_value ();
370
371 tmp = map.getfield ("v0");
372 if (tmp.is_defined ())
373 {
374 if (a_is_complex || b_is_complex)
376 else
377 resid = ColumnVector (tmp.vector_value ());
378 }
379
380 tmp = map.getfield ("disp");
381 if (tmp.is_defined ())
382 disp = tmp.strict_int_value ();
383
384 tmp = map.getfield ("cholB");
385 if (tmp.is_defined ())
386 {
387 if (tmp.numel () != 1)
388 error ("eigs: OPTS.cholB must be a scalar value");
389
390 cholB = tmp.strict_bool_value ("eigs: OPTS.cholB must be a logical value");
391 }
392
393 tmp = map.getfield ("permB");
394 if (tmp.is_defined ())
395 permB = ColumnVector (tmp.vector_value ()) - 1.0;
396 }
397
398 if (nargin > (4+arg_offset))
399 error ("eigs: incorrect number of arguments");
400
401 // Test undeclared (no issym) matrix inputs for symmetry
402 if (! sym_tested && ! have_a_fcn)
403 {
404 if (a_is_complex)
405 {
406 if (a_is_sparse)
407 symmetric = ascm.ishermitian ();
408 else
409 symmetric = acm.ishermitian ();
410 }
411 else
412 {
413 if (a_is_sparse)
414 symmetric = asmm.issymmetric ();
415 else
416 symmetric = amm.issymmetric ();
417 }
418 }
419
420 // Matrix B
421 if (have_b)
422 {
423 if (a_is_complex || b_is_complex)
424 {
425 if (b_is_sparse)
426 {
427 bscm = args(b_arg).sparse_complex_matrix_value ();
428 // Validation is 17x faster in C++ than in m-file.
429 if (bscm.any_element_is_inf_or_nan ())
430 error ("eigs: matrix A contains Inf or NaN values");
431 }
432 else
433 {
434 bcm = args(b_arg).complex_matrix_value ();
435 if (bcm.any_element_is_inf_or_nan ())
436 error ("eigs: matrix B contains Inf or NaN values");
437 }
438 }
439 else
440 {
441 if (b_is_sparse)
442 {
443 bsmm = args(b_arg).sparse_matrix_value ();
444 if (bsmm.any_element_is_inf_or_nan ())
445 error ("eigs: matrix B contains Inf or NaN values");
446 }
447 else
448 {
449 bmm = args(b_arg).matrix_value ();
450 if (bmm.any_element_is_inf_or_nan ())
451 error ("eigs: matrix B contains Inf or NaN values");
452 }
453 }
454 }
455
456 // Mode 1 for SM mode seems unstable for some reason.
457 // Use Mode 3 instead, with sigma = 0.
458 if (! have_sigma && typ == "SM")
459 have_sigma = true;
460
461 octave_idx_type nconv;
462 if (a_is_complex || b_is_complex)
463 {
465 eigs_complex_fcn = [&callback] (const ComplexColumnVector& x,
466 int& eigs_error)
467 {
468 return callback.eigs_complex_func (x, eigs_error);
469 };
470
471 ComplexMatrix eig_vec;
472 ComplexColumnVector eig_val;
473
474 if (have_a_fcn)
475 {
476 if (b_is_sparse)
478 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
479 eig_val, bscm, permB, cresid, octave_stdout, tol,
480 (nargout > 1), cholB, disp, maxit);
481 else
483 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
484 eig_val, bcm, permB, cresid, octave_stdout, tol,
485 (nargout > 1), cholB, disp, maxit);
486 }
487 else if (have_sigma)
488 {
489 if (a_is_sparse)
491 (ascm, sigma, k, p, info, eig_vec, eig_val, bscm, permB,
492 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
493 maxit);
494 else
496 (acm, sigma, k, p, info, eig_vec, eig_val, bcm, permB,
497 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
498 maxit);
499 }
500 else
501 {
502 if (a_is_sparse)
504 (ascm, typ, k, p, info, eig_vec, eig_val, bscm, permB,
505 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
506 maxit);
507 else
509 (acm, typ, k, p, info, eig_vec, eig_val, bcm, permB,
510 cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
511 maxit);
512 }
513
514 if (nargout < 2)
515 {
516 if (symmetric)
517 retval(0) = real (eig_val);
518 else
519 retval(0) = eig_val;
520 }
521 else
522 {
523 if (symmetric)
524 retval = ovl (eig_vec, DiagMatrix (real (eig_val)), double (info));
525 else
526 retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
527 }
528 }
529 else if (sigmai != 0.0)
530 {
532 eigs_complex_fcn = [&callback] (const ComplexColumnVector& x,
533 int& eigs_error)
534 {
535 return callback.eigs_complex_func (x, eigs_error);
536 };
537
538 // Promote real problem to a complex one.
539 ComplexMatrix eig_vec;
540 ComplexColumnVector eig_val;
541
542 if (have_a_fcn)
543 {
544 if (b_is_sparse)
546 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
547 eig_val, bscm, permB, cresid, octave_stdout, tol,
548 (nargout > 1), cholB, disp, maxit);
549 else
551 (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
552 eig_val, bcm, permB, cresid, octave_stdout, tol,
553 (nargout > 1), cholB, disp, maxit);
554 }
555 else
556 {
557 if (a_is_sparse)
559 (SparseComplexMatrix (asmm), sigma, k, p, info, eig_vec,
560 eig_val, SparseComplexMatrix (bsmm), permB, cresid,
561 octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
562 else
564 (ComplexMatrix (amm), sigma, k, p, info, eig_vec,
565 eig_val, ComplexMatrix (bmm), permB, cresid,
566 octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
567 }
568
569 if (nargout < 2)
570 {
571 if (symmetric)
572 retval(0) = real (eig_val);
573 else
574 retval(0) = eig_val;
575 }
576 else
577 {
578 if (symmetric)
579 retval = ovl (eig_vec, DiagMatrix (real (eig_val)), double (info));
580 else
581 retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
582 }
583 }
584 else
585 {
586 EigsFunc eigs_fcn = [&callback] (const ColumnVector& x, int& eigs_error)
587 {
588 return callback.eigs_func (x, eigs_error);
589 };
590
591 if (symmetric)
592 {
593 Matrix eig_vec;
594 ColumnVector eig_val;
595
596 if (have_a_fcn)
597 {
598 if (b_is_sparse)
600 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
601 eig_val, bsmm, permB, resid, octave_stdout, tol,
602 (nargout > 1), cholB, disp, maxit);
603 else
605 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
606 eig_val, bmm, permB, resid, octave_stdout, tol,
607 (nargout > 1), cholB, disp, maxit);
608 }
609 else if (have_sigma)
610 {
611 if (a_is_sparse)
613 (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm,
614 permB, resid, octave_stdout, tol, (nargout > 1),
615 cholB, disp, maxit);
616 else
618 (amm, sigmar, k, p, info, eig_vec, eig_val, bmm,
619 permB, resid, octave_stdout, tol, (nargout > 1),
620 cholB, disp, maxit);
621 }
622 else
623 {
624 if (a_is_sparse)
626 (asmm, typ, k, p, info, eig_vec, eig_val, bsmm,
627 permB, resid, octave_stdout, tol, (nargout > 1),
628 cholB, disp, maxit);
629 else
631 (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
632 resid, octave_stdout, tol, (nargout > 1), cholB,
633 disp, maxit);
634 }
635
636 if (nargout < 2)
637 retval(0) = eig_val;
638 else
639 retval = ovl (eig_vec, DiagMatrix (eig_val), double (info));
640 }
641 else
642 {
643 ComplexMatrix eig_vec;
644 ComplexColumnVector eig_val;
645
646 if (have_a_fcn)
647 {
648 if (b_is_sparse)
650 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
651 eig_val, bsmm, permB, resid, octave_stdout, tol,
652 (nargout > 1), cholB, disp, maxit);
653 else
655 (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
656 eig_val, bmm, permB, resid, octave_stdout, tol,
657 (nargout > 1), cholB, disp, maxit);
658 }
659 else if (have_sigma)
660 {
661 if (a_is_sparse)
663 (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm,
664 permB, resid, octave_stdout, tol, (nargout > 1),
665 cholB, disp, maxit);
666 else
668 (amm, sigmar, k, p, info, eig_vec, eig_val, bmm,
669 permB, resid, octave_stdout, tol, (nargout > 1),
670 cholB, disp, maxit);
671 }
672 else
673 {
674 if (a_is_sparse)
676 (asmm, typ, k, p, info, eig_vec, eig_val, bsmm,
677 permB, resid, octave_stdout, tol, (nargout > 1),
678 cholB, disp, maxit);
679 else
681 (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
682 resid, octave_stdout, tol, (nargout > 1), cholB,
683 disp, maxit);
684 }
685
686 if (nargout < 2)
687 retval(0) = eig_val;
688 else
689 retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
690 }
691 }
692
693 if (nconv <= 0)
694 warning_with_id ("Octave:eigs:UnconvergedEigenvalues",
695 "eigs: None of the %" OCTAVE_IDX_TYPE_FORMAT
696 " requested eigenvalues converged", k);
697 else if (nconv < k)
698 warning_with_id ("Octave:eigs:UnconvergedEigenvalues",
699 "eigs: Only %" OCTAVE_IDX_TYPE_FORMAT
700 " of the %" OCTAVE_IDX_TYPE_FORMAT
701 " requested eigenvalues converged",
702 nconv, k);
703
704 if (! fcn_name.empty ())
705 {
706 symbol_table& symtab = interp.get_symbol_table ();
707
708 symtab.clear_function (fcn_name);
709 }
710
711 return retval;
712
713#else
714
715 octave_unused_parameter (interp);
716 octave_unused_parameter (args);
717 octave_unused_parameter (nargout);
718
719 err_disabled_feature ("eigs", "ARPACK");
720
721#endif
722}
723
724/*
725## No test needed for internal helper function.
726%!assert (1)
727*/
728
729OCTAVE_END_NAMESPACE(octave)
bool ishermitian() const
Definition CMatrix.cc:175
bool any_element_is_inf_or_nan() const
Definition CNDArray.cc:271
bool issymmetric() const
Definition dMatrix.cc:137
bool any_element_is_inf_or_nan() const
Definition dNDArray.cc:324
bool any_element_is_inf_or_nan() const
Definition CSparse.cc:7506
bool ishermitian() const
Definition CSparse.cc:140
bool issymmetric() const
Definition dSparse.cc:128
bool any_element_is_inf_or_nan() const
Definition dSparse.cc:7414
octave_value getfield(const std::string &key) const
Definition oct-map.cc:183
octave_idx_type length() const
Definition ovl.h:111
int strict_int_value(bool frc_str_conv=false) const
Definition ov.h:813
bool is_defined() const
Definition ov.h:590
octave_idx_type numel() const
Definition ov.h:557
Array< Complex > complex_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
Array< double > vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
bool strict_bool_value() const
Definition ov.h:892
double double_value(bool frc_str_conv=false) const
Definition ov.h:845
void clear_function(const std::string &name)
Definition symtab.cc:446
ColumnVector real(const ComplexColumnVector &a)
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage()
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 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)
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)
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)
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)
octave_idx_type EigsRealSymmetricFunc(EigsFunc fcn, 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)
octave_idx_type EigsRealNonSymmetricFunc(EigsFunc fcn, 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)
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)
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:791
octave_idx_type EigsComplexNonSymmetricFunc(EigsComplexFunc fcn, 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)
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:1083
void warning_with_id(const char *id, const char *fmt,...)
Definition error.cc:1098
void error(const char *fmt,...)
Definition error.cc:1008
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:5709
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:217
#define octave_stdout
Definition pager.h:301
F77_RET_T const F77_DBLE * x