1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD>
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc/slepc.css" type="text/css">
<TITLE>EPSSetEigenvalueComparison</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<H1>EPSSetEigenvalueComparison</H1>
Specifies the eigenvalue comparison function when <A HREF="../EPS/EPSSetWhichEigenpairs.html#EPSSetWhichEigenpairs">EPSSetWhichEigenpairs</A>() is set to EPS_WHICH_USER.
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepceps.h"
PetscErrorCode EPSSetEigenvalueComparison(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
</PRE>
Logically Collective on <A HREF="../EPS/EPS.html#EPS">EPS</A>
<P>
<H3><FONT COLOR="#883300">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>eps </B></TD><TD> - eigensolver context obtained from <A HREF="../EPS/EPSCreate.html#EPSCreate">EPSCreate</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>func </B></TD><TD> - a pointer to the comparison function
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD> - a context pointer (the last parameter to the comparison function)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Calling Sequence of func</FONT></H3>
<pre>
func(<A HREF="../EPS/EPS.html#EPS">EPS</A> eps,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
</pre>
<P>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>eps </B></TD><TD> - eigensolver context obtained from <A HREF="../EPS/EPSCreate.html#EPSCreate">EPSCreate</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ar </B></TD><TD> - real part of the 1st eigenvalue
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ai </B></TD><TD> - imaginary part of the 1st eigenvalue
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>br </B></TD><TD> - real part of the 2nd eigenvalue
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>bi </B></TD><TD> - imaginary part of the 2nd eigenvalue
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>res </B></TD><TD> - result of comparison
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD> - optional context, as set by <A HREF="../EPS/EPSSetEigenvalueComparison.html#EPSSetEigenvalueComparison">EPSSetEigenvalueComparison</A>()
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Note</FONT></H3>
<H3><FONT COLOR="#883300">The returning parameter 'res' can be</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>negative </B></TD><TD> - if the 1st eigenvalue is preferred to the 2st one
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>zero </B></TD><TD> - if both eigenvalues are equally preferred
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>positive </B></TD><TD> - if the 2st eigenvalue is preferred to the 1st one
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
<A HREF="../EPS/EPSSetWhichEigenpairs.html#EPSSetWhichEigenpairs">EPSSetWhichEigenpairs</A>(), <A HREF="../EPS/EPSSortEigenvalues.html#EPSSortEigenvalues">EPSSortEigenvalues</A>(), <A HREF="../EPS/EPSWhich.html#EPSWhich">EPSWhich</A>
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/eps/interface/opts.c.html#EPSSetEigenvalueComparison">src/eps/interface/opts.c</A>
<BR><A HREF="./index.html">Index of all EPS routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
<P><H3><FONT COLOR="#883300">Examples</FONT></H3>
<A HREF="../../../src/eps/examples/tutorials/ex18.c.html">src/eps/examples/tutorials/ex18.c.html</A><BR>
</BODY></HTML>
|