File: NEPSetRefine.html

package info (click to toggle)
slepc 3.7.3%2Bdfsg1-5
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 23,788 kB
  • ctags: 100,221
  • sloc: ansic: 79,324; makefile: 3,897; python: 2,734; fortran: 1,139; f90: 225; sh: 100
file content (80 lines) | stat: -rw-r--r-- 4,497 bytes parent folder | download
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
77
78
79
80
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://slepc.upv.es/documentation/current/docs/manualpages/NEP/NEPSetRefine.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>NEPSetRefine</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/NEP/NEPSetRefine.html "><small>Report Typos and Errors</small></a></div>

<H1>NEPSetRefine</H1>
Specifies the refinement type (and options) to be used after the solve. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcnep.h" 
PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
</PRE>
Logically Collective on <A HREF="../NEP/NEP.html#NEP">NEP</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>nep    </B></TD><TD>&nbsp;- the nonlinear eigensolver context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>refine </B></TD><TD>&nbsp;- refinement type
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>npart  </B></TD><TD>&nbsp;- number of partitions of the communicator
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>tol    </B></TD><TD>&nbsp;- the convergence tolerance
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>its    </B></TD><TD>&nbsp;- maximum number of refinement iterations
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>scheme </B></TD><TD>&nbsp;- which scheme to be used for solving the involved linear systems
</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>-nep_refine &lt;type&gt; </B></TD><TD>&nbsp;- refinement type, one of &lt;none,simple,multiple&gt;
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-nep_refine_partitions &lt;n&gt; </B></TD><TD>&nbsp;- the number of partitions
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-nep_refine_tol &lt;tol&gt; </B></TD><TD>&nbsp;- the tolerance
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-nep_refine_its &lt;its&gt; </B></TD><TD>&nbsp;- number of iterations
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-nep_refine_scheme </B></TD><TD>&nbsp;- to set the scheme for the linear solves
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
By default, iterative refinement is disabled, since it may be very
costly. There are two possible refinement strategies: simple and multiple.
The simple approach performs iterative refinement on each of the
converged eigenpairs individually, whereas the multiple strategy works
with the invariant pair as a whole, refining all eigenpairs simultaneously.
The latter may be required for the case of multiple eigenvalues.
<P>
In some cases, especially when using direct solvers within the
iterative refinement method, it may be helpful for improved scalability
to split the communicator in several partitions. The npart parameter
indicates how many partitions to use (defaults to 1).
<P>
The tol and its parameters specify the stopping criterion. In the simple
method, refinement continues until the residual of each eigenpair is
below the tolerance (tol defaults to the <A HREF="../NEP/NEP.html#NEP">NEP</A> tol, but may be set to a
different value). In contrast, the multiple method simply performs its
refinement iterations (just one by default).
<P>
The scheme argument is used to change the way in which linear systems are
solved. Possible choices are: explicit, mixed block elimination (MBE),
and Schur complement.
<P>

<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
 <A HREF="../NEP/NEPGetRefine.html#NEPGetRefine">NEPGetRefine</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/nep/interface/nepopts.c.html#NEPSetRefine">src/nep/interface/nepopts.c</A>
<BR><A HREF="./index.html">Index of all NEP 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>