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
|
<!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>EPSSortDenseSchur</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<H1>EPSSortDenseSchur</H1>
Reorders the Schur decomposition computed by <A HREF="../EPS/EPSDenseSchur.html#EPSDenseSchur">EPSDenseSchur</A>().
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepceps.h"
PetscErrorCode EPSSortDenseSchur(EPS eps,PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
</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>eps </B></TD><TD> - the eigensolver context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>n </B></TD><TD> - dimension of the matrix
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>k </B></TD><TD> - first active column
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ldt </B></TD><TD> - leading dimension of T
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Input/Output Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>T </B></TD><TD> - the upper (quasi-)triangular matrix
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Q </B></TD><TD> - the orthogonal matrix of Schur vectors
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>wr </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></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
This function reorders the eigenvalues in wr,wi located in positions k
to n according to the sort order specified in EPSetWhicheigenpairs.
The Schur decomposition Z*T*Z^T, is also reordered by means of rotations
so that eigenvalues in the diagonal blocks of T follow the same order.
<P>
Both T and Q are overwritten.
<P>
This routine uses LAPACK routines xTREXC.
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
<A HREF="../EPS/EPSDenseHessenberg.html#EPSDenseHessenberg">EPSDenseHessenberg</A>(), <A HREF="../EPS/EPSDenseSchur.html#EPSDenseSchur">EPSDenseSchur</A>(), <A HREF="../EPS/EPSDenseTridiagonal.html#EPSDenseTridiagonal">EPSDenseTridiagonal</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/eps/interface/dense.c.html#EPSSortDenseSchur">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>
|