slepc-3.20.2 2024-03-15
Report Typos and Errors

NEPApplyFunction

Applies the nonlinear function T(lambda) to a given vector.

Synopsis

#include "slepcnep.h" 
PetscErrorCode NEPApplyFunction(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A,Mat B)
Collective

Input Parameters

nep  - the nonlinear eigensolver context
lambda  - scalar argument
x  - vector to be multiplied against
v  - workspace vector (used only in the case of split form)

Output Parameters

y  - result vector
A  - (optional) Function matrix, for callback interface only
B  - (unused) preconditioning matrix

Note

If the nonlinear operator is represented in split form, the result y = T(lambda)*x is computed without building T(lambda) explicitly. In that case, parameters A and B are not used. Otherwise, the matrix T(lambda) is built and the effect is the same as a call to NEPComputeFunction() followed by a MatMult().

See Also

NEPSetSplitOperator(), NEPComputeFunction(), NEPApplyAdjoint()

Level

developer

Location

src/nep/interface/nepsolve.c

Index of all NEP routines
Table of Contents for all manual pages
Index of all manual pages