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 70 71 72 73 74 75 76
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://slepc.upv.es/documentation/current/docs/manualpages/PEP/PEPSetScale.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>PEPSetScale</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>slepc-3.7.3 2016-09-29</b></div>
<div id="bugreport" align=right><a href="mailto:slepc-maint@upv.es?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: slepc-3.7.3 v3.7.3 docs/manualpages/PEP/PEPSetScale.html "><small>Report Typos and Errors</small></a></div>
<H1>PEPSetScale</H1>
Specifies the scaling strategy to be used.
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcpep.h"
PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
</PRE>
Logically Collective on <A HREF="../PEP/PEP.html#PEP">PEP</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>pep </B></TD><TD> - the eigensolver context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>scale </B></TD><TD> - scaling strategy
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>alpha </B></TD><TD> - the scaling factor used in the scalar strategy
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Dl </B></TD><TD> - the left diagonal matrix of the diagonal scaling algorithm
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Dr </B></TD><TD> - the right diagonal matrix of the diagonal scaling algorithm
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>its </B></TD><TD> - number of iterations of the diagonal scaling algorithm
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>lambda </B></TD><TD> - approximation to wanted eigenvalues (modulus)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pep_scale <type> </B></TD><TD> - scaling type, one of <none,scalar,diagonal,both>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pep_scale_factor <alpha> </B></TD><TD> - the scaling factor
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pep_scale_its <its> </B></TD><TD> - number of iterations
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pep_scale_lambda <lambda> </B></TD><TD> - approximation to eigenvalues
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
There are two non-exclusive scaling strategies: scalar and diagonal.
<P>
In the scalar strategy, scaling is applied to the eigenvalue, that is,
mu = lambda/alpha is the new eigenvalue and all matrices are scaled
accordingly. After solving the scaled problem, the original lambda is
recovered. Parameter 'alpha' must be positive. Use PETSC_DECIDE to let
the solver compute a reasonable scaling factor.
<P>
In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
of the computed results in some cases. The user may provide the Dr and Dl
matrices represented as Vec objects storing diagonal elements. If not
provided, these matrices are computed internally. This option requires
that the polynomial coefficient matrices are of MATAIJ type.
The parameter 'its' is the number of iterations performed by the method.
Parameter 'lambda' must be positive. Use PETSC_DECIDE or set lambda = 1.0 if
no information about eigenvalues is available.
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
<A HREF="../PEP/PEPGetScale.html#PEPGetScale">PEPGetScale</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/pep/interface/pepopts.c.html#PEPSetScale">src/pep/interface/pepopts.c</A>
<BR><A HREF="./index.html">Index of all PEP 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>
|