slepc-main 2024-11-15
Report Typos and Errors

BVTraceQuadrature

Computes an estimate of the number of eigenvalues inside a region via quantities computed in the quadrature rule of contour integral methods.

Synopsis

#include "slepcbv.h" 
PetscErrorCode BVTraceQuadrature(BV Y,BV V,PetscInt L,PetscInt L_max,PetscScalar *w,VecScatter scat,PetscSubcomm subcomm,PetscInt npoints,PetscBool useconj,PetscReal *est_eig)
Collective

Input Parameters

Y  - first basis vectors
V  - second basis vectors
L  - block size
L_max  - maximum block size
w  - quadrature weights
scat  - (optional) VecScatter object to communicate between subcommunicators
subcomm  - subcommunicator layout
npoints  - number of points to process by the subcommunicator
useconj  - whether conjugate points can be used or not

Output Parameter

est_eig  - estimated eigenvalue count

Notes

This function returns an estimation of the number of eigenvalues in the region, computed as trace(V'*S_0), where S_0 is the first panel of S computed by BVSumQuadrature().

When using subcommunicators, Y is stored in the subcommunicators for a subset of integration points. In that case, the computation is done in the subcomm and then scattered to the whole communicator in S using the VecScatter scat. The value npoints is the number of points to be processed in this subcomm and the flag useconj indicates whether symmetric points can be reused.

See Also

BVScatter(), BVDotQuadrature(), BVSumQuadrature(), RGComputeQuadrature(), RGCanUseConjugates()

Level

developer

Location

src/sys/classes/bv/interface/bvcontour.c

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