slepc-main 2024-12-17
Report Typos and Errors

RGCheckInside

Determines if a set of given points are inside the region or not.

Synopsis

#include "slepcrg.h" 
PetscErrorCode RGCheckInside(RG rg,PetscInt n,PetscScalar *ar,PetscScalar *ai,PetscInt *inside)
Not Collective

Input Parameters

rg  - the region context
n  - number of points to check
ar  - array of real parts
ai  - array of imaginary parts

Output Parameter

inside  - array of results (1=inside, 0=on the contour, -1=outside)

Note

The point a is expressed as a couple of PetscScalar variables ar,ai. If built with complex scalars, the point is supposed to be stored in ar, otherwise ar,ai contain the real and imaginary parts, respectively.

If a scaling factor was set, the points are scaled before checking.

See Also

RGSetScale(), RGSetComplement()

Level

intermediate

Location

src/sys/classes/rg/interface/rgbasic.c

Index of all RG routines
Table of Contents for all manual pages
Index of all manual pages