EVSL  1.1.0
EigenValues Slicing Library
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros
Functions
ratfilter.c File Reference

Computing and applying rational filters. More...

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <complex.h>
#include <math.h>
#include "def.h"
#include "blaslapack.h"
#include "struct.h"
#include "internal_proto.h"

Go to the source code of this file.

Functions

void contQuad (int method, int n, complex double *zk)
 Compute the locations of the poles. More...
 
void ratf2p2 (int n, int *mulp, complex double *zk, complex double *alp, int m, double *z, double *xx)
 Compute the function value of the multiple pole rational filter at real locations. More...
 
void pfe2 (complex double s1, complex double s2, int k1, int k2, complex double *alp, complex double *bet)
 Get the fraction expansion of 1/[(z-s1)^k1 (z-s2)^k2]. More...
 
complex double integg2 (complex double s1, complex double s2, complex double *alp, int k1, complex double *bet, int k2, double a, double b)
 Integration of 1/[(z-s1)^k1 (z-s2)^k2] from a to b. More...
 
void weights (int n, complex double *zk, int *mulp, double lambda, complex double *omega)
 Compute the LS weight for each multiple pole. More...
 
int scaleweigthts (int n, double a, double b, complex double *zk, int *mulp, complex double *omegaM)
 Compute the weights and pole locations on [a, b]. More...
 
void set_ratf_def (ratparams *rat)
 Sets default values for ratparams struct. More...
 
int find_ratf (double *intv, ratparams *rat)
 
void free_rat (ratparams *rat)
 
void RatFiltApply (int n, ratparams *rat, double *b, double *x, double *w6)
 Apply rational filter R to a vetor b. More...
 

Detailed Description

Computing and applying rational filters.

Definition in file ratfilter.c.

Function Documentation

void contQuad ( int  method,
int  n,
complex double *  zk 
)

Compute the locations of the poles.

----------------------— Cauchy integration-based filter -----------—

Parameters
method0 for Guass Legendre; 1 for midpoint
nNumber of poles in the upper half plane
[out]zkVector of pole locations

Definition at line 23 of file ratfilter.c.

References DSTEV(), M_PI, and Malloc.

Referenced by find_ratf().

Here is the call graph for this function:

Here is the caller graph for this function:

int find_ratf ( double *  intv,
ratparams rat 
)

Parameters
intv= an array of length 4 [intv[0], intv[1]] is the interval of desired eigenvalues [intv[2], intv[3]] is the global interval of all eigenvalues it must contain all eigenvalues of A

OUT:

Parameters
[out]ratthese are set in rat struct:
omega : expansion coefficients of rational filter
zk : location of the poles used
aa : adjusted left endpoint
bb : adjusted right endpoint
dd : half-width and..
cc : ..center of interval

Definition at line 375 of file ratfilter.c.

References _ratparams::aa, _ratparams::ASIGBsol, _ratparams::bb, _ratparams::beta, contQuad(), Malloc, max, _ratparams::method, min, _ratparams::mulp, _ratparams::num, _ratparams::omega, _ratparams::pow, _ratparams::pw, scaleweigthts(), weights(), and _ratparams::zk.

Referenced by evsl_find_rat(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

void free_rat ( ratparams rat)

Definition at line 424 of file ratfilter.c.

References _ratparams::ASIGBsol, _ratparams::mulp, _ratparams::omega, and _ratparams::zk.

Referenced by evsl_free_rat(), and main().

Here is the caller graph for this function:

complex double integg2 ( complex double  s1,
complex double  s2,
complex double *  alp,
int  k1,
complex double *  bet,
int  k2,
double  a,
double  b 
)

Integration of 1/[(z-s1)^k1 (z-s2)^k2] from a to b.

Definition at line 129 of file ratfilter.c.

Referenced by weights().

Here is the caller graph for this function:

void pfe2 ( complex double  s1,
complex double  s2,
int  k1,
int  k2,
complex double *  alp,
complex double *  bet 
)

Get the fraction expansion of 1/[(z-s1)^k1 (z-s2)^k2].

Definition at line 95 of file ratfilter.c.

Referenced by weights().

Here is the caller graph for this function:

void ratf2p2 ( int  n,
int *  mulp,
complex double *  zk,
complex double *  alp,
int  m,
double *  z,
double *  xx 
)

Compute the function value of the multiple pole rational filter at real locations.

---------------—Multiple pole rational filter evaluation -----------—

Parameters
nnumber of the pole
mulpmultiplicity of the pole
zkarray containing the poles.
alpfractional expansion coefficients
mnumber of locations to be evaluated
zreal locations to be evaluated
[out]xx: function values at real locations z

Definition at line 69 of file ratfilter.c.

Referenced by scaleweigthts(), and weights().

Here is the caller graph for this function:

void RatFiltApply ( int  n,
ratparams rat,
double *  b,
double *  x,
double *  w6 
)

Apply rational filter R to a vetor b.

For generalized e.v problem x = L' * (A-SB) \ L*b [w:=work] x = L * b w = (A-sB) \ x x = L' * w

Parameters
[in]ratratparams struct
[in]nLength of array
[in]bx = L * b
w6Work array of size 4*n for standard ev problem, size 6*n for generalized ev problem
[out]xBecomes R(A)b

Definition at line 448 of file ratfilter.c.

References _ratparams::ASIGBsol, DAXPY(), DSCAL(), evsl_timer(), evsldata, evslstat, _evsldata::ifGenEv, _ratparams::mulp, _evslstat::n_ratAv, _ratparams::num, _ratparams::omega, and _evslstat::t_ratAv.

Referenced by RatLanNr(), and RatLanTr().

Here is the call graph for this function:

Here is the caller graph for this function:

int scaleweigthts ( int  n,
double  a,
double  b,
complex double *  zk,
int *  mulp,
complex double *  omegaM 
)

Compute the weights and pole locations on [a, b].

---------------—Transform poles and weights computed on [-1, 1] to [a, b] -------—

Parameters
nnumber of poles used in the upper half plane
a,b[a, b] is the interval of desired eigenvalues
zklocation of the poles
mulpmultiplicity of the poles
[out]omegaMmultiple LS weights

Definition at line 313 of file ratfilter.c.

References ratf2p2().

Referenced by find_ratf().

Here is the call graph for this function:

Here is the caller graph for this function:

void set_ratf_def ( ratparams rat)

Sets default values for ratparams struct.

Definition at line 345 of file ratfilter.c.

References _ratparams::aa, _ratparams::bar, _ratparams::bb, _ratparams::beta, _ratparams::method, _ratparams::num, and _ratparams::pw.

Referenced by evsl_find_rat(), and main().

Here is the caller graph for this function:

void weights ( int  n,
complex double *  zk,
int *  mulp,
double  lambda,
complex double *  omega 
)

Compute the LS weight for each multiple pole.

---------------—multiple pole LS rational filter weights-----------—

Parameters
nnumber of poles in the upper half plane
zkpole locations
mulpmultiplicity of each pole
lambdaLS integration weight for [-1, 1]
[out]omegaLS weight for each pole

Definition at line 169 of file ratfilter.c.

References integg2(), Malloc, pfe2(), ratf2p2(), and ZGESV().

Referenced by find_ratf().

Here is the call graph for this function:

Here is the caller graph for this function: