GNU Octave  9.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-2024 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 
52 struct eigs_callback
53 {
54 public:
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?
78 static int call_depth = 0;
79 
81 eigs_callback::eigs_func (const ColumnVector& x, int& eigs_error)
82 {
83  ColumnVector retval;
84  octave_value_list args;
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 
121 eigs_callback::eigs_complex_func (const ComplexColumnVector& x,
122  int& eigs_error)
123 {
124  ComplexColumnVector retval;
125  octave_value_list args;
126  args(0) = x;
127 
128  if (m_eigs_fcn.is_defined ())
129  {
130  octave_value_list tmp;
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 
157 DEFMETHOD (__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{})
179 Undocumented 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  callback.m_eigs_fcn = get_function_handle (interp, args(0), "x");
234 
235  if (callback.m_eigs_fcn.is_undefined ())
236  error ("eigs: unknown function");
237 
238  if (nargin < 2)
239  error ("eigs: incorrect number of arguments");
240 
241  n = args(1).nint_value ();
242  arg_offset = 1;
243  have_a_fcn = true;
244  }
245  else
246  {
247  if (args(0).iscomplex ())
248  {
249  if (args(0).issparse ())
250  {
251  ascm = (args(0).sparse_complex_matrix_value ());
252  a_is_sparse = true;
253  }
254  else
255  acm = (args(0).complex_matrix_value ());
256  a_is_complex = true;
257  }
258  else
259  {
260  if (args(0).issparse ())
261  {
262  asmm = (args(0).sparse_matrix_value ());
263  a_is_sparse = true;
264  }
265  else
266  {
267  amm = (args(0).matrix_value ());
268  }
269  }
270  }
271 
272  // Note hold off reading B until later to avoid issues of double
273  // copies of the matrix if B is full/real while A is complex.
274  if (nargin > 1 + arg_offset
275  && ! (args(1 + arg_offset).is_real_scalar ()))
276  {
277  if (args(1+arg_offset).iscomplex ())
278  {
279  b_arg = 1+arg_offset;
280  if (args(b_arg).issparse ())
281  {
282  bscm = (args(b_arg).sparse_complex_matrix_value ());
283  b_is_sparse = true;
284  }
285  else
286  bcm = (args(b_arg).complex_matrix_value ());
287  have_b = true;
288  b_is_complex = true;
289  arg_offset++;
290  }
291  else
292  {
293  b_arg = 1+arg_offset;
294  if (args(b_arg).issparse ())
295  {
296  bsmm = (args(b_arg).sparse_matrix_value ());
297  b_is_sparse = true;
298  }
299  else
300  bmm = (args(b_arg).matrix_value ());
301  have_b = true;
302  arg_offset++;
303  }
304  }
305 
306  if (nargin > (1+arg_offset))
307  k = args(1+arg_offset).nint_value ();
308 
309  if (nargin > (2+arg_offset))
310  {
311  if (args(2+arg_offset).is_string ())
312  {
313  typ = args(2+arg_offset).string_value ();
314 
315  // Use STL function to convert to upper case
316  transform (typ.begin (), typ.end (), typ.begin (), toupper);
317 
318  sigma = 0.0;
319  }
320  else
321  {
322  sigma = args(2+arg_offset).xcomplex_value ("eigs: SIGMA must be a scalar or a string");
323 
324  have_sigma = true;
325  }
326  }
327 
328  sigmar = sigma.real ();
329  sigmai = sigma.imag ();
330 
331  if (nargin > (3+arg_offset))
332  {
333  if (! args(3+arg_offset).isstruct ())
334  error ("eigs: OPTS argument must be a structure");
335 
336  octave_scalar_map map = args(3
337  +arg_offset).xscalar_map_value ("eigs: OPTS argument must be a scalar structure");
338 
339  octave_value tmp;
340 
341  // issym is ignored for complex matrix inputs
342  tmp = map.getfield ("issym");
343  if (tmp.is_defined ())
344  {
345  if (tmp.numel () != 1)
346  error ("eigs: OPTS.issym must be a scalar value");
347 
348  symmetric = tmp.xbool_value ("eigs: OPTS.issym must be a logical value");
349  sym_tested = true;
350  }
351 
352  // isreal is ignored if A is not a function
353  if (have_a_fcn)
354  {
355  tmp = map.getfield ("isreal");
356  if (tmp.is_defined ())
357  {
358  if (tmp.numel () != 1)
359  error ("eigs: OPTS.isreal must be a scalar value");
360 
361  a_is_complex = ! tmp.xbool_value ("eigs: OPTS.isreal must be a logical value");
362  }
363  }
364 
365  tmp = map.getfield ("tol");
366  if (tmp.is_defined ())
367  tol = tmp.double_value ();
368 
369  tmp = map.getfield ("maxit");
370  if (tmp.is_defined ())
371  maxit = tmp.nint_value ();
372 
373  tmp = map.getfield ("p");
374  if (tmp.is_defined ())
375  p = tmp.nint_value ();
376 
377  tmp = map.getfield ("v0");
378  if (tmp.is_defined ())
379  {
380  if (a_is_complex || b_is_complex)
381  cresid = ComplexColumnVector (tmp.complex_vector_value ());
382  else
383  resid = ColumnVector (tmp.vector_value ());
384  }
385 
386  tmp = map.getfield ("disp");
387  if (tmp.is_defined ())
388  disp = tmp.nint_value ();
389 
390  tmp = map.getfield ("cholB");
391  if (tmp.is_defined ())
392  {
393  if (tmp.numel () != 1)
394  error ("eigs: OPTS.cholB must be a scalar value");
395 
396  cholB = tmp.xbool_value ("eigs: OPTS.cholB must be a logical value");
397  }
398 
399  tmp = map.getfield ("permB");
400  if (tmp.is_defined ())
401  permB = ColumnVector (tmp.vector_value ()) - 1.0;
402  }
403 
404  if (nargin > (4+arg_offset))
405  error ("eigs: incorrect number of arguments");
406 
407  // Test undeclared (no issym) matrix inputs for symmetry
408  if (! sym_tested && ! have_a_fcn)
409  {
410  if (a_is_complex)
411  {
412  if (a_is_sparse)
413  symmetric = ascm.ishermitian ();
414  else
415  symmetric = acm.ishermitian ();
416  }
417  else
418  {
419  if (a_is_sparse)
420  symmetric = asmm.issymmetric ();
421  else
422  symmetric = amm.issymmetric ();
423  }
424  }
425 
426  if (have_b)
427  {
428  if (a_is_complex || b_is_complex)
429  {
430  if (b_is_sparse)
431  bscm = args(b_arg).sparse_complex_matrix_value ();
432  else
433  bcm = args(b_arg).complex_matrix_value ();
434  }
435  else
436  {
437  if (b_is_sparse)
438  bsmm = args(b_arg).sparse_matrix_value ();
439  else
440  bmm = args(b_arg).matrix_value ();
441  }
442  }
443 
444  // Mode 1 for SM mode seems unstable for some reason.
445  // Use Mode 3 instead, with sigma = 0.
446  if (! have_sigma && typ == "SM")
447  have_sigma = true;
448 
449  octave_idx_type nconv;
450  if (a_is_complex || b_is_complex)
451  {
453  eigs_complex_fcn = [&callback] (const ComplexColumnVector& x,
454  int& eigs_error)
455  {
456  return callback.eigs_complex_func (x, eigs_error);
457  };
458 
459  ComplexMatrix eig_vec;
460  ComplexColumnVector eig_val;
461 
462  if (have_a_fcn)
463  {
464  if (b_is_sparse)
466  (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
467  eig_val, bscm, permB, cresid, octave_stdout, tol,
468  (nargout > 1), cholB, disp, maxit);
469  else
471  (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
472  eig_val, bcm, permB, cresid, octave_stdout, tol,
473  (nargout > 1), cholB, disp, maxit);
474  }
475  else if (have_sigma)
476  {
477  if (a_is_sparse)
479  (ascm, sigma, k, p, info, eig_vec, eig_val, bscm, permB,
480  cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
481  maxit);
482  else
484  (acm, sigma, k, p, info, eig_vec, eig_val, bcm, permB,
485  cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
486  maxit);
487  }
488  else
489  {
490  if (a_is_sparse)
492  (ascm, typ, k, p, info, eig_vec, eig_val, bscm, permB,
493  cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
494  maxit);
495  else
497  (acm, typ, k, p, info, eig_vec, eig_val, bcm, permB,
498  cresid, octave_stdout, tol, (nargout > 1), cholB, disp,
499  maxit);
500  }
501 
502  if (nargout < 2)
503  {
504  if (symmetric)
505  retval(0) = real (eig_val);
506  else
507  retval(0) = eig_val;
508  }
509  else
510  {
511  if (symmetric)
512  retval = ovl (eig_vec, DiagMatrix (real (eig_val)), double (info));
513  else
514  retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
515  }
516  }
517  else if (sigmai != 0.0)
518  {
520  eigs_complex_fcn = [&callback] (const ComplexColumnVector& x,
521  int& eigs_error)
522  {
523  return callback.eigs_complex_func (x, eigs_error);
524  };
525 
526  // Promote real problem to a complex one.
527  ComplexMatrix eig_vec;
528  ComplexColumnVector eig_val;
529 
530  if (have_a_fcn)
531  {
532  if (b_is_sparse)
534  (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
535  eig_val, bscm, permB, cresid, octave_stdout, tol,
536  (nargout > 1), cholB, disp, maxit);
537  else
539  (eigs_complex_fcn, n, typ, sigma, k, p, info, eig_vec,
540  eig_val, bcm, permB, cresid, octave_stdout, tol,
541  (nargout > 1), cholB, disp, maxit);
542  }
543  else
544  {
545  if (a_is_sparse)
547  (SparseComplexMatrix (asmm), sigma, k, p, info, eig_vec,
548  eig_val, SparseComplexMatrix (bsmm), permB, cresid,
549  octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
550  else
552  (ComplexMatrix (amm), sigma, k, p, info, eig_vec,
553  eig_val, ComplexMatrix (bmm), permB, cresid,
554  octave_stdout, tol, (nargout > 1), cholB, disp, maxit);
555  }
556 
557  if (nargout < 2)
558  {
559  if (symmetric)
560  retval(0) = real (eig_val);
561  else
562  retval(0) = eig_val;
563  }
564  else
565  {
566  if (symmetric)
567  retval = ovl (eig_vec, DiagMatrix (real (eig_val)), double (info));
568  else
569  retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
570  }
571  }
572  else
573  {
574  EigsFunc eigs_fcn = [&callback] (const ColumnVector& x, int& eigs_error)
575  {
576  return callback.eigs_func (x, eigs_error);
577  };
578 
579  if (symmetric)
580  {
581  Matrix eig_vec;
582  ColumnVector eig_val;
583 
584  if (have_a_fcn)
585  {
586  if (b_is_sparse)
587  nconv = EigsRealSymmetricFunc
588  (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
589  eig_val, bsmm, permB, resid, octave_stdout, tol,
590  (nargout > 1), cholB, disp, maxit);
591  else
592  nconv = EigsRealSymmetricFunc
593  (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
594  eig_val, bmm, permB, resid, octave_stdout, tol,
595  (nargout > 1), cholB, disp, maxit);
596  }
597  else if (have_sigma)
598  {
599  if (a_is_sparse)
601  (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm,
602  permB, resid, octave_stdout, tol, (nargout > 1),
603  cholB, disp, maxit);
604  else
606  (amm, sigmar, k, p, info, eig_vec, eig_val, bmm,
607  permB, resid, octave_stdout, tol, (nargout > 1),
608  cholB, disp, maxit);
609  }
610  else
611  {
612  if (a_is_sparse)
614  (asmm, typ, k, p, info, eig_vec, eig_val, bsmm,
615  permB, resid, octave_stdout, tol, (nargout > 1),
616  cholB, disp, maxit);
617  else
619  (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
620  resid, octave_stdout, tol, (nargout > 1), cholB,
621  disp, maxit);
622  }
623 
624  if (nargout < 2)
625  retval(0) = eig_val;
626  else
627  retval = ovl (eig_vec, DiagMatrix (eig_val), double (info));
628  }
629  else
630  {
631  ComplexMatrix eig_vec;
632  ComplexColumnVector eig_val;
633 
634  if (have_a_fcn)
635  {
636  if (b_is_sparse)
638  (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
639  eig_val, bsmm, permB, resid, octave_stdout, tol,
640  (nargout > 1), cholB, disp, maxit);
641  else
643  (eigs_fcn, n, typ, sigmar, k, p, info, eig_vec,
644  eig_val, bmm, permB, resid, octave_stdout, tol,
645  (nargout > 1), cholB, disp, maxit);
646  }
647  else if (have_sigma)
648  {
649  if (a_is_sparse)
651  (asmm, sigmar, k, p, info, eig_vec, eig_val, bsmm,
652  permB, resid, octave_stdout, tol, (nargout > 1),
653  cholB, disp, maxit);
654  else
656  (amm, sigmar, k, p, info, eig_vec, eig_val, bmm,
657  permB, resid, octave_stdout, tol, (nargout > 1),
658  cholB, disp, maxit);
659  }
660  else
661  {
662  if (a_is_sparse)
664  (asmm, typ, k, p, info, eig_vec, eig_val, bsmm,
665  permB, resid, octave_stdout, tol, (nargout > 1),
666  cholB, disp, maxit);
667  else
669  (amm, typ, k, p, info, eig_vec, eig_val, bmm, permB,
670  resid, octave_stdout, tol, (nargout > 1), cholB,
671  disp, maxit);
672  }
673 
674  if (nargout < 2)
675  retval(0) = eig_val;
676  else
677  retval = ovl (eig_vec, ComplexDiagMatrix (eig_val), double (info));
678  }
679  }
680 
681  if (nconv <= 0)
682  warning_with_id ("Octave:eigs:UnconvergedEigenvalues",
683  "eigs: None of the %" OCTAVE_IDX_TYPE_FORMAT
684  " requested eigenvalues converged", k);
685  else if (nconv < k)
686  warning_with_id ("Octave:eigs:UnconvergedEigenvalues",
687  "eigs: Only %" OCTAVE_IDX_TYPE_FORMAT
688  " of the %" OCTAVE_IDX_TYPE_FORMAT
689  " requested eigenvalues converged",
690  nconv, k);
691 
692  if (! fcn_name.empty ())
693  {
694  symbol_table& symtab = interp.get_symbol_table ();
695 
696  symtab.clear_function (fcn_name);
697  }
698 
699  return retval;
700 
701 #else
702 
703  octave_unused_parameter (interp);
704  octave_unused_parameter (args);
705  octave_unused_parameter (nargout);
706 
707  err_disabled_feature ("eigs", "ARPACK");
708 
709 #endif
710 }
711 
712 /*
713 ## No test needed for internal helper function.
714 %!assert (1)
715 */
716 
717 OCTAVE_END_NAMESPACE(octave)
bool ishermitian() const
Definition: CMatrix.cc:174
Definition: dMatrix.h:42
bool issymmetric() const
Definition: dMatrix.cc:136
bool ishermitian() const
Definition: CSparse.cc:142
bool issymmetric() const
Definition: dSparse.cc:131
octave_value getfield(const std::string &key) const
Definition: oct-map.cc:183
octave_idx_type length() const
Definition: ovl.h:113
bool xbool_value(const char *fmt,...) const
bool is_defined() const
Definition: ov.h:592
octave_idx_type numel() const
Definition: ov.h:559
int nint_value(bool frc_str_conv=false) const
Definition: ov.h:819
Array< double > vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
Array< Complex > complex_vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
double double_value(bool frc_str_conv=false) const
Definition: ov.h:841
void clear_function(const std::string &name)
Definition: symtab.cc:446
ColumnVector real(const ComplexColumnVector &a)
Definition: dColVector.cc:137
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
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 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 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)
Definition: eigs-base.cc:1381
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)
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
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)
Definition: eigs-base.cc:3553
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:1063
void warning_with_id(const char *id, const char *fmt,...)
Definition: error.cc:1078
void() error(const char *fmt,...)
Definition: error.cc:988
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:5468
octave_value get_function_handle(interpreter &interp, const octave_value &arg, const std::string &parameter_name)
F77_RET_T const F77_DBLE * x
octave_idx_type n
Definition: mx-inlines.cc:761
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:219
#define octave_stdout
Definition: pager.h:309