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
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="https://slepc.upv.es/documentation/current//Users/jroman/tmp/slepc-3.23.1/docs/manualpages/PEP/PEPSTOARGetInertias.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>PEPSTOARGetInertias</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>slepc-3.23.1 2025-05-01</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.23.1 v3.23.1 /Users/jroman/tmp/slepc-3.23.1/docs/manualpages/PEP/PEPSTOARGetInertias.html "><small>Report Typos and Errors</small></a></div>
<H1>PEPSTOARGetInertias</H1>
Gets the values of the shifts and their corresponding inertias in case of doing spectrum slicing for a computational interval.
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcpep.h"
<A HREF="https://petsc.org/release/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../PEP/PEPSTOARGetInertias.html#PEPSTOARGetInertias">PEPSTOARGetInertias</A>(<A HREF="../PEP/PEP.html#PEP">PEP</A> pep,<A HREF="https://petsc.org/release/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> *n,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> *shifts[],<A HREF="https://petsc.org/release/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> *inertias[]) PeNS
</PRE>
Not Collective
<P>
<H3><FONT COLOR="#883300">Input Parameter</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 eigenproblem solver context
</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>n </B></TD><TD> - number of shifts, including the endpoints of the interval
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>shifts </B></TD><TD> - the values of the shifts used internally in the solver
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>inertias </B></TD><TD> - the values of the inertia in each shift
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
If called after <A HREF="../PEP/PEPSolve.html#PEPSolve">PEPSolve</A>(), all shifts used internally by the solver are
returned (including both endpoints and any intermediate ones). If called
before <A HREF="../PEP/PEPSolve.html#PEPSolve">PEPSolve</A>() and after <A HREF="../PEP/PEPSetUp.html#PEPSetUp">PEPSetUp</A>() then only the information of the
endpoints of subintervals is available.
<P>
This function is only available for spectrum slicing runs.
<P>
The returned arrays should be freed by the user. Can pass NULL in any of
the two arrays if not required.
<P>
<H3><FONT COLOR="#883300">Fortran Notes</FONT></H3>
The calling sequence from Fortran is
<PRE>
<A HREF="../PEP/PEPSTOARGetInertias.html#PEPSTOARGetInertias">PEPSTOARGetInertias</A>(pep,n,shifts,inertias,ierr)
integer n
double precision shifts(*)
integer inertias(*)
</PRE>
The arrays should be at least of length n. The value of n can be determined
by an initial call
<PRE>
<A HREF="../PEP/PEPSTOARGetInertias.html#PEPSTOARGetInertias">PEPSTOARGetInertias</A>(pep,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
</PRE>
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
<A HREF="../PEP/PEPSetInterval.html#PEPSetInterval">PEPSetInterval</A>()
<BR><P><B></B><H3><FONT COLOR="#883300">Level</FONT></H3>advanced<BR>
<H3><FONT COLOR="#883300">Location</FONT></H3>
</B><A HREF="../../../src/pep/impls/krylov/stoar/stoar.c.html#PEPSTOARGetInertias">src/pep/impls/krylov/stoar/stoar.c</A>
<P><H3><FONT COLOR="#883300">Examples</FONT></H3>
<A HREF="../../../src/pep/tutorials/ex38.c.html">src/pep/tutorials/ex38.c</A><BR>
<BR><BR><A HREF="./index.html">Index of all PEP routines</A>
<BR><A HREF="../../../docs/manual.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|