![]() |
EVSL
1.1.0
EigenValues Slicing Library
|
Default solver function for solving shifted systems and factoring B. More...
#include <string.h>#include "def.h"#include "struct.h"#include "internal_proto.h"#include "cholmod.h"#include "umfpack.h"#include "evsl_direct.h"Go to the source code of this file.
Data Structures | |
| struct | _BSolDataDirect |
Macros | |
| #define | CHOLMOD(name) cholmod_ ## name |
Typedefs | |
| typedef struct _BSolDataDirect | BSolDataDirect |
Functions | |
| cholmod_sparse * | csrMat_to_cholmod_sparse (csrMat *A, int stype) |
| Create cholmod_sparse matrix just as a wrapper of a csrMat This will be useful for the B matrix, since it will be factored. More... | |
| void | arr_to_cholmod_dense (int m, int n, double *A, cholmod_dense *B) |
| wrap an array in cholmod dense (no copying) More... | |
| void | BSolDirect (double *b, double *x, void *data) |
| Solver function of B. More... | |
| int | SetupBSolDirect (csrMat *B, void **data) |
| Setup the B-sol by computing the Cholesky factorization of B. More... | |
| void | FreeBSolDirectData (void *vdata) |
| Free solver data. More... | |
| void | LTSolDirect (double *b, double *x, void *data) |
| Solver function of L^{T} x = L^{-T}*b. More... | |
| void | ASIGMABSolDirect (int n, double *br, double *bi, double *xr, double *xz, void *data) |
| default complex linear solver routine passed to evsl More... | |
| int | SetupASIGMABSolDirect (csrMat *A, csrMat *BB, int num, complex double *zk, void **data) |
| setup the default solver for A - SIGMA B More... | |
| void | FreeASIGMABSolDirect (int num, void **data) |
| free the data needed by UMFpack More... | |
Variables | |
| cholmod_dense | cholmod_X |
| cholmod_dense | cholmod_B |
| cholmod_dense * | cholmod_Y =NULL |
| cholmod_dense * | cholmod_E =NULL |
| cholmod_dense * | cholmod_W =NULL |
Default solver function for solving shifted systems and factoring B.
This file contains the default solver function for solving shifted systems with A - SIGMA*I or A - SIGMA*B, and for factoring the B matrix. The default solvers are UMFPACK and CHOLMOD from SuiteSparse
Definition in file evsl_suitesparse.c.
| #define CHOLMOD | ( | name | ) | cholmod_ ## name |
Definition at line 13 of file evsl_suitesparse.c.
Referenced by BSolDirect(), FreeBSolDirectData(), LTSolDirect(), and SetupBSolDirect().
| typedef struct _BSolDataDirect BSolDataDirect |
| void arr_to_cholmod_dense | ( | int | m, |
| int | n, | ||
| double * | A, | ||
| cholmod_dense * | B | ||
| ) |
wrap an array in cholmod dense (no copying)
Definition at line 79 of file evsl_suitesparse.c.
Referenced by SetupBSolDirect().

| void ASIGMABSolDirect | ( | int | n, |
| double * | br, | ||
| double * | bi, | ||
| double * | xr, | ||
| double * | xz, | ||
| void * | data | ||
| ) |
default complex linear solver routine passed to evsl
complex linear solver routine passed to evsl
| n | size of the system | |
| br,bi | vectors of length n, complex right-hand side (real and imaginary) | |
| data | all data that are needed for solving the system | |
| [out] | xr,xz | vectors of length n, complex solution (real and imaginary) |
Definition at line 249 of file evsl_suitesparse.c.
References _ASBSolDataDirect::b, CHKERR, cs_ci_ipvec(), cs_ci_lsolve(), cs_ci_usolve(), cs_complex_t, cs_ci_numeric::L, _ASBSolDataDirect::n, _ASBSolDataDirect::N, cs_ci_numeric::pinv, cs_ci_symbolic::q, _ASBSolDataDirect::S, cs_ci_numeric::U, and _ASBSolDataDirect::x.

| void BSolDirect | ( | double * | b, |
| double * | x, | ||
| void * | data | ||
| ) |
Solver function of B.
Definition at line 118 of file evsl_suitesparse.c.
References CHOLMOD, cholmod_B, cholmod_E, cholmod_X, cholmod_Y, _BSolDataDirect::cm, cs_di_ipvec(), cs_di_lsolve(), cs_di_ltsolve(), cs_di_pvec(), cs_di_numeric::L, _BSolDataDirect::LB, _BSolDataDirect::n, _BSolDataDirect::N, cs_di_symbolic::pinv, _BSolDataDirect::S, and _BSolDataDirect::w.

| cholmod_sparse* csrMat_to_cholmod_sparse | ( | csrMat * | A, |
| int | stype | ||
| ) |
Create cholmod_sparse matrix just as a wrapper of a csrMat This will be useful for the B matrix, since it will be factored.
Definition at line 34 of file evsl_suitesparse.c.
References _csrMat::a, _csrMat::ia, _csrMat::ja, Malloc, _csrMat::ncols, _csrMat::nrows, and SuiteSparse_long.
Referenced by SetupBSolDirect().

| void FreeASIGMABSolDirect | ( | int | num, |
| void ** | data | ||
| ) |
free the data needed by UMFpack
free the data needed by CXSparse
Definition at line 380 of file evsl_suitesparse.c.
References _ASBSolDataDirect::b, cs_ci_nfree(), cs_ci_sfree(), _ASBSolDataDirect::N, _ASBSolDataDirect::S, and _ASBSolDataDirect::x.

| void FreeBSolDirectData | ( | void * | vdata | ) |
Free solver data.
Definition at line 198 of file evsl_suitesparse.c.
References CHOLMOD, cholmod_E, cholmod_W, cholmod_Y, _BSolDataDirect::cm, cs_di_nfree(), cs_di_sfree(), _BSolDataDirect::LB, _BSolDataDirect::N, _BSolDataDirect::S, and _BSolDataDirect::w.

| void LTSolDirect | ( | double * | b, |
| double * | x, | ||
| void * | data | ||
| ) |
Solver function of L^{T} x = L^{-T}*b.
Definition at line 215 of file evsl_suitesparse.c.
References CHOLMOD, cholmod_B, cholmod_E, cholmod_W, cholmod_X, cholmod_Y, _BSolDataDirect::cm, cs_di_ltsolve(), cs_di_pvec(), cs_di_numeric::L, _BSolDataDirect::LB, _BSolDataDirect::n, _BSolDataDirect::N, cs_di_symbolic::pinv, _BSolDataDirect::S, and _BSolDataDirect::w.

setup the default solver for A - SIGMA B
setup CXsparse solver for A - SIGMA B
The setup invovles shifting the matrix A - SIGMA B and factorizing the shifted matrix The solver function and the data will be saved data Generally speaking, each pole can have a different solver
| A | matrix A |
| BB | matrix B, if NULL, it means B is identity |
| num | the number of SIGMA's |
| zk | array of SIGMA's of length num |
| data | all data that are needed for solving the system |
Definition at line 272 of file evsl_suitesparse.c.
References _csrMat::a, _ASBSolDataDirect::b, Calloc, CHKERR, cs_ci_lu(), cs_ci_sqr(), cs_complex_t, evsl_timer(), evslstat, free_csr(), cs_ci_sparse::i, _csrMat::ia, _csrMat::ja, cs_ci_sparse::m, Malloc, matadd(), _ASBSolDataDirect::n, _ASBSolDataDirect::N, cs_ci_sparse::n, _csrMat::ncols, _csrMat::nrows, cs_ci_sparse::nz, cs_ci_sparse::nzmax, cs_ci_sparse::p, _ASBSolDataDirect::S, speye(), SuiteSparse_long, _evslstat::t_setASigBsv, _ASBSolDataDirect::x, and cs_ci_sparse::x.

| int SetupBSolDirect | ( | csrMat * | B, |
| void ** | data | ||
| ) |
Setup the B-sol by computing the Cholesky factorization of B.
| B | matrix B |
| data | Struct which will be initialized |
Definition at line 144 of file evsl_suitesparse.c.
References _csrMat::a, arr_to_cholmod_dense(), CHKERR, CHOLMOD, cholmod_B, cholmod_X, _BSolDataDirect::cm, cs_di_chol(), cs_di_schol(), cs_di_spfree(), csrMat_to_cholmod_sparse(), evsl_timer(), evslstat, cs_di_sparse::i, _csrMat::ia, _csrMat::ja, _BSolDataDirect::LB, cs_di_sparse::m, Malloc, _BSolDataDirect::n, _BSolDataDirect::N, cs_di_sparse::n, _csrMat::nrows, cs_di_sparse::nz, cs_di_sparse::nzmax, cs_di_sparse::p, _BSolDataDirect::S, _evslstat::t_setBsv, _BSolDataDirect::w, and cs_di_sparse::x.

| cholmod_dense cholmod_B |
Definition at line 112 of file evsl_suitesparse.c.
Referenced by BSolDirect(), LTSolDirect(), and SetupBSolDirect().
| cholmod_dense * cholmod_E =NULL |
Definition at line 112 of file evsl_suitesparse.c.
Referenced by BSolDirect(), FreeBSolDirectData(), and LTSolDirect().
| cholmod_dense * cholmod_W =NULL |
Definition at line 113 of file evsl_suitesparse.c.
Referenced by FreeBSolDirectData(), and LTSolDirect().
| cholmod_dense cholmod_X |
global variables for cholmod solve
Definition at line 112 of file evsl_suitesparse.c.
Referenced by BSolDirect(), LTSolDirect(), and SetupBSolDirect().
| cholmod_dense * cholmod_Y =NULL |
Definition at line 112 of file evsl_suitesparse.c.
Referenced by BSolDirect(), FreeBSolDirectData(), and LTSolDirect().
1.8.6