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 81 82 83 84 85 86
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="https://slepc.upv.es/documentation/current/docs/manualpages/RG/RGRingSetParameters.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>RGRingSetParameters</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>slepc-3.22.2 2024-12-02</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.22.2 v3.22.2 docs/manualpages/RG/RGRingSetParameters.html "><small>Report Typos and Errors</small></a></div>
<H1>RGRingSetParameters</H1>
Sets the parameters defining the ring region.
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcrg.h"
<A HREF="https://petsc.org/release/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../RG/RGRingSetParameters.html#RGRingSetParameters">RGRingSetParameters</A>(<A HREF="../RG/RG.html#RG">RG</A> rg,<A HREF="https://petsc.org/release/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> center,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> radius,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> vscale,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> start_ang,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> end_ang,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> width)
</PRE>
Logically 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>rg </B></TD><TD> - the region context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>center </B></TD><TD> - center of the ellipse
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>radius </B></TD><TD> - radius of the ellipse
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>vscale </B></TD><TD> - vertical scale of the ellipse
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>start_ang </B></TD><TD> - the right-hand side angle
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>end_ang </B></TD><TD> - the left-hand side angle
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>width </B></TD><TD> - width of the ring
</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>-rg_ring_center </B></TD><TD> - Sets the center
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-rg_ring_radius </B></TD><TD> - Sets the radius
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-rg_ring_vscale </B></TD><TD> - Sets the vertical scale
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-rg_ring_startangle </B></TD><TD> - Sets the right-hand side angle
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-rg_ring_endangle </B></TD><TD> - Sets the left-hand side angle
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-rg_ring_width </B></TD><TD> - Sets the width of the ring
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
The values of center, radius and vscale have the same meaning as in the
ellipse region. The startangle and endangle define the span of the ring
(by default it is the whole ring), while the width is the separation
between the two concentric ellipses (above and below the radius by
width/2).
<P>
The start and end angles are expressed as a fraction of the circumference.
The allowed range is [0..1], with 0 corresponding to 0 radians, 0.25 to
pi/2 radians, and so on. It is allowed to have startangle>endangle, in
which case the ring region crosses over the zero angle.
<P>
In the case of complex scalars, a complex center can be provided in the
command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
-rg_ring_center 1.0+2.0i
<P>
When PETSc is built with real scalars, the center is restricted to a real value,
and the start and end angles must be such that the region is symmetric with
respect to the real axis.
<P>
For all arguments except center, you can use <A HREF="https://petsc.org/release/manualpages/Sys/PETSC_CURRENT.html#PETSC_CURRENT">PETSC_CURRENT</A> to keep the current
value, and <A HREF="https://petsc.org/release/manualpages/Sys/PETSC_DETERMINE.html#PETSC_DETERMINE">PETSC_DETERMINE</A> to set them to a default value (1 for radius, vscale,
end_ang, 0 for start_and, and 0.1 for width).
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
<A HREF="../RG/RGRingGetParameters.html#RGRingGetParameters">RGRingGetParameters</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/sys/classes/rg/impls/ring/rgring.c.html#RGRingSetParameters">src/sys/classes/rg/impls/ring/rgring.c</A>
<BR><BR><A HREF="./index.html">Index of all RG 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>
|