| 12
 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
 
 | <!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>EPSDenseGNHEP</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<H1>EPSDenseGNHEP</H1>
Solves a dense Generalized non-Hermitian Eigenvalue Problem. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepceps.h" 
PetscErrorCode EPSDenseGNHEP(PetscInt n_,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
</PRE>
Not Collective
<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>n  </B></TD><TD> - dimension of the eigenproblem
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>A  </B></TD><TD> - pointer to the array containing the matrix values for A
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>B  </B></TD><TD> - pointer to the array containing the matrix values for B
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Output Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>w  </B></TD><TD> - pointer to the array to store the computed eigenvalues
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>wi </B></TD><TD> - imaginary part of the eigenvalues (only when using real numbers)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>V  </B></TD><TD> - pointer to the array to store right eigenvectors
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>W  </B></TD><TD> - pointer to the array to store left eigenvectors
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
If either V or W are PETSC_NULL then the corresponding eigenvectors are
not computed.
<P>
Matrices A and B are overwritten.
<P>
This routine uses LAPACK routines xGGEVX.
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
 <A HREF="../EPS/EPSDenseNHEP.html#EPSDenseNHEP">EPSDenseNHEP</A>(), <A HREF="../EPS/EPSDenseHEP.html#EPSDenseHEP">EPSDenseHEP</A>(), <A HREF="../EPS/EPSDenseGHEP.html#EPSDenseGHEP">EPSDenseGHEP</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/eps/interface/dense.c.html#EPSDenseGNHEP">src/eps/interface/dense.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>
</BODY></HTML>
 |