HLIBpro
1.2
|
Coordinates | |
hlib_coord_t | hlib_coord_import (const size_t n, const unsigned int dim, double **coord, const double *period, int *info) |
size_t | hlib_coord_ncoord (const hlib_coord_t coord, int *info) |
unsigned int | hlib_coord_dim (const hlib_coord_t coord, int *info) |
double * | hlib_coord_get (const hlib_coord_t coord, const size_t i, int *info) |
double * | hlib_coord_bbmin (const hlib_coord_t coord, const size_t i, int *info) |
double * | hlib_coord_bbmax (const hlib_coord_t coord, const size_t i, int *info) |
int | hlib_coord_has_bbox (const hlib_coord_t coord, int *info) |
void | hlib_coord_set_period (hlib_coord_t coord, const double period[], int *info) |
void | hlib_coord_get_period (hlib_coord_t coord, double period[], int *info) |
void | hlib_coord_free (hlib_coord_t coord, int *info) |
size_t | hlib_coord_bytesize (const hlib_coord_t coord, int *info) |
void | hlib_coord_print_vrml (const hlib_coord_t coord, const char *filename, int *info) |
Admissibility condition | |
hlib_admcond_t | hlib_admcond_geom (const hlib_adm_t crit, const double eta, int *info) |
hlib_admcond_t | hlib_admcond_geom_period (const hlib_adm_t crit, const double eta, const double *period, const unsigned int dim, int *info) |
hlib_admcond_t | hlib_admcond_alg (const hlib_adm_t crit, const double eta, const hlib_matrix_t S, const hlib_permutation_t row_perm, const hlib_permutation_t col_perm, int *info) |
void | hlib_admcond_free (const hlib_admcond_t ac, int *info) |
Clusters | |
int | hlib_cl_first (const hlib_cluster_t cl, int *info) |
int | hlib_cl_last (const hlib_cluster_t cl, int *info) |
size_t | hlib_cl_size (const hlib_cluster_t cl, int *info) |
size_t | hlib_cl_nsons (const hlib_cluster_t cl, int *info) |
hlib_cluster_t | hlib_cl_son (const hlib_cluster_t cl, const unsigned int i, int *info) |
int | hlib_cl_is_leaf (const hlib_cluster_t cl, int *info) |
size_t | hlib_cl_nnodes (const hlib_cluster_t cl, int *info) |
size_t | hlib_cl_depth (const hlib_cluster_t cl, int *info) |
hlib_cluster_t | hlib_cl_copy (const hlib_cluster_t cl, int *info) |
void | hlib_cl_free (hlib_cluster_t cl, int *info) |
size_t | hlib_cl_bytesize (const hlib_cluster_t cl, int *info) |
void | hlib_cl_print_ps (const hlib_cluster_t cl, const char *filename, int *info) |
Cluster Trees | |
hlib_clustertree_t | hlib_clt_build_bsp (const hlib_coord_t coord, const hlib_bsp_t bsptype, const unsigned int nmin, int *info) |
hlib_clustertree_t | hlib_clt_build_bsp_nd (const hlib_coord_t coord, const hlib_matrix_t S, const hlib_bsp_t bsptype, const unsigned int nmin, int *info) |
hlib_clustertree_t | hlib_clt_build_bsp_part (const hlib_coord_t coord, const unsigned int *partition, const int eq_depth, const hlib_bsp_t bsptype, const unsigned int nmin, int *info) |
hlib_clustertree_t | hlib_clt_build_alg (const hlib_matrix_t S, const unsigned int nmin, int *info) |
hlib_clustertree_t | hlib_clt_build_alg_nd (const hlib_matrix_t S, const unsigned int nmin, int *info) |
hlib_clustertree_t | hlib_clt_build_alg_part (const hlib_matrix_t S, const unsigned int *partition, const unsigned int nmin, int *info) |
hlib_cluster_t | hlib_clt_root (hlib_clustertree_t ct, int *info) |
hlib_permutation_t | hlib_clt_perm_i2e (hlib_clustertree_t ct, int *info) |
hlib_permutation_t | hlib_clt_perm_e2i (hlib_clustertree_t ct, int *info) |
void | hlib_clt_free (hlib_clustertree_t ct, int *info) |
size_t | hlib_clt_bytesize (const hlib_clustertree_t ct, int *info) |
size_t | hlib_clt_nnodes (const hlib_clustertree_t ct, int *info) |
size_t | hlib_clt_depth (const hlib_clustertree_t ct, int *info) |
void | hlib_clt_print_ps (const hlib_clustertree_t ct, const char *filename, int *info) |
Block Clusters | |
hlib_blockcluster_t | hlib_bc_parent (const hlib_blockcluster_t bc, int *info) |
size_t | hlib_bc_nsons (const hlib_blockcluster_t bc, int *info) |
hlib_blockcluster_t | hlib_bc_son (const hlib_blockcluster_t bc, const unsigned int i, const unsigned int j, int *info) |
int | hlib_bc_is_leaf (const hlib_blockcluster_t bc, int *info) |
int | hlib_bc_is_adm (const hlib_blockcluster_t bc, int *info) |
size_t | hlib_bc_nnodes (const hlib_blockcluster_t bc, int *info) |
size_t | hlib_bc_depth (const hlib_blockcluster_t bc, int *info) |
hlib_cluster_t | hlib_bc_rowcl (const hlib_blockcluster_t bc, int *info) |
hlib_cluster_t | hlib_bc_colcl (const hlib_blockcluster_t bc, int *info) |
hlib_blockcluster_t | hlib_bc_copy (const hlib_blockcluster_t bc, int *info) |
void | hlib_bc_free (hlib_blockcluster_t bc, int *info) |
size_t | hlib_bc_bytesize (const hlib_blockcluster_t bc, int *info) |
int | hlib_bc_csp (const hlib_blockcluster_t bc, int *info) |
void | hlib_bc_print_ps (const hlib_blockcluster_t bc, const char *filename, int *info) |
Block Cluster Trees | |
Root of block cluster trees with access to index permutations. | |
hlib_blockclustertree_t | hlib_bct_build (const hlib_clustertree_t rowct, const hlib_clustertree_t colct, const hlib_admcond_t ac, int *info) |
size_t | hlib_bct_nnodes (const hlib_blockclustertree_t bct, int *info) |
size_t | hlib_bct_depth (const hlib_blockclustertree_t bct, int *info) |
hlib_clustertree_t | hlib_bct_rowct (const hlib_blockclustertree_t bct, int *info) |
hlib_clustertree_t | hlib_bct_colct (const hlib_blockclustertree_t bct, int *info) |
void | hlib_bct_free (hlib_blockclustertree_t bct, int *info) |
size_t | hlib_bct_bytesize (const hlib_blockclustertree_t bct, int *info) |
void | hlib_bct_distribute_block (hlib_blockclustertree_t bct, const unsigned int p, const unsigned int min_lvl, const int symmetric, int *info) |
int | hlib_bct_csp (const hlib_blockclustertree_t bct, int *info) |
void | hlib_bct_print_ps (const hlib_blockclustertree_t bct, const char *filename, int *info) |
Vectors | |
size_t | hlib_vector_size (const hlib_vector_t v, int *info) |
hlib_vector_t | hlib_vector_copy (const hlib_vector_t v, int *info) |
size_t | hlib_vector_bytesize (const hlib_vector_t v, int *info) |
void | hlib_vector_free (hlib_vector_t v, int *info) |
hlib_real_t | hlib_vector_entry_get (const hlib_vector_t x, const size_t i, int *info) |
hlib_complex_t | hlib_vector_centry_get (const hlib_vector_t x, const size_t i, int *info) |
void | hlib_vector_entry_set (const hlib_vector_t x, const size_t i, const hlib_real_t f, int *info) |
void | hlib_vector_centry_set (const hlib_vector_t x, const size_t i, const hlib_complex_t f, int *info) |
hlib_vector_t | hlib_vector_build (const size_t size, int *info) |
hlib_vector_t | hlib_vector_cbuild (const size_t size, int *info) |
void | hlib_vector_import (hlib_vector_t v, const hlib_real_t *arr, int *info) |
void | hlib_vector_cimport (hlib_vector_t v, const hlib_complex_t *arr, int *info) |
void | hlib_vector_export (const hlib_vector_t v, hlib_real_t *arr, int *info) |
void | hlib_vector_cexport (const hlib_vector_t v, hlib_complex_t *arr, int *info) |
void | hlib_vector_fill (hlib_vector_t x, const hlib_real_t f, int *info) |
void | hlib_vector_cfill (hlib_vector_t x, const hlib_complex_t f, int *info) |
void | hlib_vector_fill_rand (hlib_vector_t x, int *info) |
void | hlib_vector_assign (hlib_vector_t y, const hlib_vector_t x, int *info) |
void | hlib_vector_scale (hlib_vector_t x, const hlib_real_t f, int *info) |
void | hlib_vector_cscale (hlib_vector_t x, const hlib_complex_t f, int *info) |
void | hlib_vector_axpy (const hlib_real_t alpha, const hlib_vector_t x, hlib_vector_t y, int *info) |
void | hlib_vector_caxpy (const hlib_complex_t alpha, const hlib_vector_t x, hlib_vector_t y, int *info) |
hlib_complex_t | hlib_vector_dot (const hlib_vector_t x, const hlib_vector_t y, int *info) |
hlib_real_t | hlib_vector_norm2 (const hlib_vector_t x, int *info) |
hlib_real_t | hlib_vector_norm_inf (const hlib_vector_t x, int *info) |
hlib_vector_t | hlib_vector_restrict_re (const hlib_vector_t v, int *info) |
hlib_vector_t | hlib_vector_restrict_im (const hlib_vector_t v, int *info) |
void | hlib_vector_fft (const hlib_vector_t v, int *info) |
void | hlib_vector_ifft (const hlib_vector_t v, int *info) |
Matrices | |
size_t | hlib_matrix_rows (const hlib_matrix_t A, int *info) |
size_t | hlib_matrix_cols (const hlib_matrix_t A, int *info) |
hlib_matrix_t | hlib_matrix_copy (const hlib_matrix_t A, int *info) |
hlib_matrix_t | hlib_matrix_copy_acc (const hlib_matrix_t A, const hlib_acc_t acc, const int coarsen, int *info) |
void | hlib_matrix_copyto (const hlib_matrix_t A, hlib_matrix_t B, int *info) |
void | hlib_matrix_copyto_acc (const hlib_matrix_t A, hlib_matrix_t B, const hlib_acc_t acc, const int coarsen, int *info) |
hlib_matrix_t | hlib_matrix_copy_blockdiag (const hlib_matrix_t A, const unsigned int lvl, int *info) |
hlib_matrix_t | hlib_matrix_copy_blockdiag_acc (const hlib_matrix_t A, const unsigned int lvl, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_copy_nearfield (const hlib_matrix_t A, const int without_rk, int *info) |
size_t | hlib_matrix_bytesize (const hlib_matrix_t A, int *info) |
void | hlib_matrix_free (hlib_matrix_t A, int *info) |
int | hlib_matrix_is_sym (const hlib_matrix_t A, int *info) |
int | hlib_matrix_is_herm (const hlib_matrix_t A, int *info) |
int | hlib_matrix_is_complex (const hlib_matrix_t A, int *info) |
void | hlib_matrix_to_real (const hlib_matrix_t A, const int force, int *info) |
void | hlib_matrix_to_complex (const hlib_matrix_t A, const int force, int *info) |
hlib_real_t | hlib_matrix_entry_get (const hlib_matrix_t A, const size_t i, const size_t j, int *info) |
hlib_complex_t | hlib_matrix_centry_get (const hlib_matrix_t A, const size_t i, const size_t j, int *info) |
hlib_vector_t | hlib_matrix_row_vector (const hlib_matrix_t A, int *info) |
hlib_vector_t | hlib_matrix_col_vector (const hlib_matrix_t A, int *info) |
hlib_matrix_t | hlib_matrix_import_crs (const size_t rows, const size_t cols, const size_t nnz, const int *rowind, const int *colind, const hlib_real_t *coeffs, const int sym, int *info) |
hlib_matrix_t | hlib_matrix_import_ccrs (const size_t rows, const size_t cols, const size_t nnz, const int *rowind, const int *colind, const hlib_complex_t *coeffs, const int sym, int *info) |
hlib_matrix_t | hlib_matrix_import_dense (const size_t rows, const size_t cols, const hlib_real_t *D, const int sym, int *info) |
hlib_matrix_t | hlib_matrix_import_cdense (const size_t rows, const size_t cols, const hlib_complex_t *D, const int sym, int *info) |
hlib_matrix_t | hlib_matrix_build_sparse (const hlib_blockclustertree_t bct, const hlib_matrix_t S, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_build_dense (const hlib_blockclustertree_t bct, const hlib_matrix_t D, const hlib_lrapx_t lrapx, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_build_coeff (const hlib_blockclustertree_t bct, const hlib_coeff_t f, void *arg, const hlib_lrapx_t lrapx, const hlib_acc_t acc, const int sym, int *info) |
hlib_matrix_t | hlib_matrix_build_ccoeff (const hlib_blockclustertree_t bct, const hlib_ccoeff_t f, void *arg, const hlib_lrapx_t lrapx, const hlib_acc_t acc, const int sym, int *info) |
hlib_matrix_t | hlib_matrix_build_identity (const hlib_blockclustertree_t bct, int *info) |
hlib_matrix_t | hlib_matrix_assemble_block (const size_t brows, const size_t bcols, hlib_matrix_t *submat, const int copy, int *info) |
void | hlib_matrix_print_ps (const hlib_matrix_t A, const char *filename, const int options, int *info) |
Algebra | |
void | hlib_matrix_mulvec (const hlib_real_t alpha, const hlib_matrix_t A, const hlib_vector_t x, const hlib_real_t beta, hlib_vector_t y, const hlib_matop_t matop, int *info) |
void | hlib_matrix_cmulvec (const hlib_complex_t alpha, const hlib_matrix_t A, const hlib_vector_t x, const hlib_complex_t beta, hlib_vector_t y, const hlib_matop_t matop, int *info) |
void | hlib_matrix_transpose (const hlib_matrix_t A, int *info) |
void | hlib_matrix_conjugate (const hlib_matrix_t A, int *info) |
void | hlib_matrix_scale (const hlib_real_t f, const hlib_matrix_t A, int *info) |
void | hlib_matrix_cscale (const hlib_complex_t f, const hlib_matrix_t A, int *info) |
void | hlib_matrix_add (const hlib_real_t alpha, const hlib_matrix_t A, const hlib_real_t beta, hlib_matrix_t B, const hlib_acc_t acc, int *info) |
void | hlib_matrix_cadd (const hlib_complex_t alpha, const hlib_matrix_t A, const hlib_complex_t beta, hlib_matrix_t B, const hlib_acc_t acc, int *info) |
void | hlib_matrix_mul (const hlib_real_t alpha, const hlib_matop_t matop_A, const hlib_matrix_t A, const hlib_matop_t matop_B, const hlib_matrix_t B, const hlib_real_t beta, hlib_matrix_t C, const hlib_acc_t acc, int *info) |
void | hlib_matrix_cmul (const hlib_complex_t alpha, const hlib_matop_t matop_A, const hlib_matrix_t A, const hlib_matop_t matop_B, const hlib_matrix_t B, const hlib_complex_t beta, hlib_matrix_t C, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_inv (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
hlib_vector_t | hlib_matrix_inv_diag (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_factorise (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_factorise_inv (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_lu (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_ldl (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
hlib_matrix_t | hlib_matrix_ll (hlib_matrix_t A, const hlib_acc_t acc, int *info) |
void | hlib_matrix_solve_lower_left (const hlib_matrix_t L, hlib_matrix_t A, const hlib_acc_t acc, const int pointwise, const int unit_diag, int *info) |
void | hlib_matrix_solve_lower_right (hlib_matrix_t A, const hlib_matrix_t L, const hlib_matop_t op_L, const hlib_acc_t acc, const int pointwise, const int unit_diag, int *info) |
void | hlib_matrix_solve_upper_right (hlib_matrix_t A, const hlib_matrix_t U, const hlib_acc_t acc, const int pointwise, const int unit_diag, int *info) |
void | hlib_matrix_solve_diag_left (const hlib_matrix_t D, const hlib_matop_t op_D, hlib_matrix_t A, const hlib_acc_t acc, const int pointwise, int *info) |
void | hlib_matrix_solve_diag_right (hlib_matrix_t A, const hlib_matrix_t D, const hlib_matop_t op_D, const hlib_acc_t acc, const int pointwise, int *info) |
hlib_real_t | hlib_matrix_norm_frobenius (const hlib_matrix_t A, int *info) |
hlib_real_t | hlib_matrix_norm_frobenius_diff (const hlib_matrix_t A, const hlib_matrix_t B, int *info) |
hlib_real_t | hlib_matrix_norm_spectral (const hlib_matrix_t A, int *info) |
hlib_real_t | hlib_matrix_norm_spectral_inv (const hlib_matrix_t A, int *info) |
hlib_real_t | hlib_matrix_norm_spectral_diff (const hlib_matrix_t A, const hlib_matrix_t B, int *info) |
hlib_real_t | hlib_matrix_norm_inv_approx (const hlib_matrix_t A, const hlib_matrix_t B, int *info) |
Solver | |
hlib_solver_t | hlib_solver_auto (int *info) |
hlib_solver_t | hlib_solver_richardson (int *info) |
hlib_solver_t | hlib_solver_cg (int *info) |
hlib_solver_t | hlib_solver_bicgstab (int *info) |
hlib_solver_t | hlib_solver_minres (int *info) |
hlib_solver_t | hlib_solver_gmres (const int restart, int *info) |
void | hlib_solver_solve (const hlib_solver_t solver, const hlib_matrix_t A, hlib_vector_t x, const hlib_vector_t b, const hlib_matrix_t W, hlib_solve_info_t *solve_info, int *info) |
void | hlib_solver_stopcrit (hlib_solver_t solver, const int maxit, const hlib_real_t abs_red, const hlib_real_t rel_red, int *info) |
void | hlib_solver_free (hlib_solver_t solver, int *info) |
Input/Output | |
hlib_matrix_t | hlib_hformat_load_matrix (const char *filename, int *info) |
hlib_vector_t | hlib_hformat_load_vector (const char *filename, int *info) |
void | hlib_hformat_save_matrix (const hlib_matrix_t A, const char *filename, int *info) |
void | hlib_hformat_save_vector (const hlib_vector_t v, const char *filename, int *info) |
hlib_coord_t | hlib_hformat_load_coord (const char *filename, int *info) |
void | hlib_hformat_save_coord (const hlib_coord_t coord, const char *filename, int *info) |
hlib_matrix_t | hlib_samg_load_matrix (const char *filename, int *info) |
void | hlib_samg_save_matrix (const hlib_matrix_t A, const char *filename, int *info) |
hlib_vector_t | hlib_samg_load_vector (const char *filename, int *info) |
void | hlib_samg_save_vector (const hlib_vector_t x, const char *filename, int *info) |
void | hlib_samg_save_coord (const hlib_coord_t coord, const char *filename, int *info) |
hlib_coord_t | hlib_samg_load_coord (const char *filename, int *info) |
hlib_matrix_t | hlib_matlab_load_matrix (const char *filename, const char *matname, int *info) |
void | hlib_matlab_save_matrix (const hlib_matrix_t M, const char *filename, const char *matname, int *info) |
hlib_vector_t | hlib_matlab_load_vector (const char *filename, const char *vecname, int *info) |
void | hlib_matlab_save_vector (const hlib_vector_t v, const char *filename, const char *vecname, int *info) |
hlib_matrix_t | hlib_pltmg_load_matrix (const char *filename, int *info) |
hlib_matrix_t | hlib_hb_load_matrix (const char *filename, int *info) |
void | hlib_hb_save_matrix (const hlib_matrix_t M, const char *filename, int *info) |
hlib_matrix_t | hlib_mm_load_matrix (const char *filename, int *info) |
hlib_matrix_t | hlib_load_matrix (const char *filename, int *info) |
hlib_vector_t | hlib_load_vector (const char *filename, int *info) |
hlib_coord_t | hlib_load_coord (const char *filename, int *info) |
Accuracy Management | |
hlib_acc_t | hlib_acc_fixed_eps (const hlib_real_t eps) |
hlib_acc_t | hlib_acc_fixed_rank (const unsigned int k) |
hlib_acc_t | hlib_acc_blocked (const hlib_acc_t *blockacc) |
Misc. | |
double | hlib_walltime () |
double | hlib_cputime () |
void | hlib_set_n_min (const unsigned int n) |
void | hlib_set_verbosity (const unsigned int verb) |
void | hlib_set_abs_eps (const hlib_real_t eps) |
void | hlib_set_coarsening (const int build, const int arith) |
void | hlib_set_recompress (const int recompress) |
void | hlib_set_diag_scale (const int scale) |
void | hlib_set_nthreads (const unsigned int p) |
void | hlib_set_progress_cb (hlib_progressfn_t fn, void *arg) |
void | hlib_set_error_fn (const hlib_errorfn_t errorfn) |
void | hlib_set_warning_fn (const hlib_errorfn_t warnfn) |
unsigned int | hlib_major_version () |
unsigned int | hlib_minor_version () |
Functions for hlib_complex_t | |
hlib_complex_t | hlib_complex (const double re, const double im) |
hlib_complex_t | hlib_cmplx_add (const hlib_complex_t a, const hlib_complex_t b) |
hlib_complex_t | hlib_cmplx_sub (const hlib_complex_t a, const hlib_complex_t b) |
hlib_complex_t | hlib_cmplx_mul (const hlib_complex_t a, const hlib_complex_t b) |
hlib_complex_t | hlib_cmplx_div (const hlib_complex_t a, const hlib_complex_t b) |
hlib_complex_t | hlib_cmplx_conj (const hlib_complex_t a) |
hlib_real_t | hlib_cmplx_abs (const hlib_complex_t a) |
hlib_complex_t | hlib_cmplx_sqrt (const hlib_complex_t a) |
hlib_complex_t | hlib_cmplx_exp (const hlib_complex_t a) |
This modules defines functions for the usage of 𝓗𝖫𝖨𝖡𝗉𝗋𝗈 in C programs instead of C++.
To include all standard 𝓗-matrix related functions add
to your source files. For using BEM functions in C, you also have to add
hlib_acc_t hlib_acc_blocked | ( | const hlib_acc_t * | blockacc | ) |
return accuracy object with blockwise accuracy defined by array blockacc of dimension blockrows × blockcolumns corresponding to first partitioning level; the array must be stored column wise
hlib_acc_t hlib_acc_fixed_eps | ( | const hlib_real_t | eps | ) |
return accuracy object with fixed accuracy eps
hlib_acc_t hlib_acc_fixed_rank | ( | const unsigned int | k | ) |
return accuracy object with fixed rank k
hlib_admcond_t hlib_admcond_alg | ( | const hlib_adm_t | crit, |
const double | eta, | ||
const hlib_matrix_t | S, | ||
const hlib_permutation_t | row_perm, | ||
const hlib_permutation_t | col_perm, | ||
int * | info | ||
) |
create admissibility condition based on algebraic connectivity in S
void hlib_admcond_free | ( | const hlib_admcond_t | ac, |
int * | info | ||
) |
free admissibility condition
hlib_admcond_t hlib_admcond_geom | ( | const hlib_adm_t | crit, |
const double | eta, | ||
int * | info | ||
) |
create admissibility condition based on geometrical data
hlib_admcond_t hlib_admcond_geom_period | ( | const hlib_adm_t | crit, |
const double | eta, | ||
const double * | period, | ||
const unsigned int | dim, | ||
int * | info | ||
) |
create admissibility condition based on geometrical data with additional periodicity in geometry defined by vector period
size_t hlib_bc_bytesize | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return size of memory in bytes used by block cluster
hlib_cluster_t hlib_bc_colcl | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return column cluster of given block cluster
hlib_blockcluster_t hlib_bc_copy | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return copy of sub tree defined by bc
int hlib_bc_csp | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
compute sparsity constant of block cluster tree
size_t hlib_bc_depth | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return depth of block cluster
void hlib_bc_free | ( | hlib_blockcluster_t | bc, |
int * | info | ||
) |
free resources coupled with block cluster tree bct
int hlib_bc_is_adm | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return 1 if block cluster is admissible and 0 otherwise
int hlib_bc_is_leaf | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return 1 if block cluster is leaf and 0 otherwise
size_t hlib_bc_nnodes | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return number of nodes in block cluster
size_t hlib_bc_nsons | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return total number of sons of block cluster bc
hlib_blockcluster_t hlib_bc_parent | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return parent block cluster of block cluster bc
void hlib_bc_print_ps | ( | const hlib_blockcluster_t | bc, |
const char * | filename, | ||
int * | info | ||
) |
print block cluster tree in PostScript format to file filename
hlib_cluster_t hlib_bc_rowcl | ( | const hlib_blockcluster_t | bc, |
int * | info | ||
) |
return row cluster of given block cluster
hlib_blockcluster_t hlib_bc_son | ( | const hlib_blockcluster_t | bc, |
const unsigned int | i, | ||
const unsigned int | j, | ||
int * | info | ||
) |
return son (i,j) of block cluster bc
hlib_blockclustertree_t hlib_bct_build | ( | const hlib_clustertree_t | rowct, |
const hlib_clustertree_t | colct, | ||
const hlib_admcond_t | ac, | ||
int * | info | ||
) |
build block cluster tree over given row and column clusters using admissibility condition ac
size_t hlib_bct_bytesize | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
return size of memory in bytes used by block cluster
hlib_clustertree_t hlib_bct_colct | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
return column cluster tree of given block cluster tree
int hlib_bct_csp | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
compute sparsity constant of block cluster tree
size_t hlib_bct_depth | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
return depth of block cluster tree
void hlib_bct_distribute_block | ( | hlib_blockclustertree_t | bct, |
const unsigned int | p, | ||
const unsigned int | min_lvl, | ||
const int | symmetric, | ||
int * | info | ||
) |
distribute block cluster tree onto p processors based on block-wise partition with unique processor per block
void hlib_bct_free | ( | hlib_blockclustertree_t | bct, |
int * | info | ||
) |
free resources coupled with cluster tree ct
size_t hlib_bct_nnodes | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
return number of nodes in block cluster tree
void hlib_bct_print_ps | ( | const hlib_blockclustertree_t | bct, |
const char * | filename, | ||
int * | info | ||
) |
print block cluster tree in PostScript format to file filename
hlib_clustertree_t hlib_bct_rowct | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
return row cluster tree of given block cluster tree
size_t hlib_cl_bytesize | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return size of memory in bytes used by cluster
hlib_cluster_t hlib_cl_copy | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return copy of sub tree defined by cl
size_t hlib_cl_depth | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return depth of cluster tree
int hlib_cl_first | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return first index in index set
void hlib_cl_free | ( | hlib_cluster_t | cl, |
int * | info | ||
) |
free resources coupled with cluster cl
int hlib_cl_is_leaf | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return 1 if cluster is leaf and 0 otherwise
int hlib_cl_last | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return last index in index set
size_t hlib_cl_nnodes | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return number of nodes in cluster tree
size_t hlib_cl_nsons | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return number of sons of cluster
void hlib_cl_print_ps | ( | const hlib_cluster_t | cl, |
const char * | filename, | ||
int * | info | ||
) |
print cluster tree in PostScript format to file filename
size_t hlib_cl_size | ( | const hlib_cluster_t | cl, |
int * | info | ||
) |
return size of index set
hlib_cluster_t hlib_cl_son | ( | const hlib_cluster_t | cl, |
const unsigned int | i, | ||
int * | info | ||
) |
return i'th son of cluster cl
hlib_clustertree_t hlib_clt_build_alg | ( | const hlib_matrix_t | S, |
const unsigned int | nmin, | ||
int * | info | ||
) |
build clustertree via algebraic clustering based on given sparse matrix S
hlib_clustertree_t hlib_clt_build_alg_nd | ( | const hlib_matrix_t | S, |
const unsigned int | nmin, | ||
int * | info | ||
) |
build clustertree via algebraic clustering with nested dissection based on given sparse matrix S
hlib_clustertree_t hlib_clt_build_alg_part | ( | const hlib_matrix_t | S, |
const unsigned int * | partition, | ||
const unsigned int | nmin, | ||
int * | info | ||
) |
build clustertree via algebraic clustering and use a predefined partition for first step
hlib_clustertree_t hlib_clt_build_bsp | ( | const hlib_coord_t | coord, |
const hlib_bsp_t | bsptype, | ||
const unsigned int | nmin, | ||
int * | info | ||
) |
build cluster tree using binary space partitioning
hlib_clustertree_t hlib_clt_build_bsp_nd | ( | const hlib_coord_t | coord, |
const hlib_matrix_t | S, | ||
const hlib_bsp_t | bsptype, | ||
const unsigned int | nmin, | ||
int * | info | ||
) |
build cluster tree using binary space partitioning and nested dissection based on given sparse matrix S
hlib_clustertree_t hlib_clt_build_bsp_part | ( | const hlib_coord_t | coord, |
const unsigned int * | partition, | ||
const int | eq_depth, | ||
const hlib_bsp_t | bsptype, | ||
const unsigned int | nmin, | ||
int * | info | ||
) |
build cluster tree using binary space partitioning but use a predefined partition for first step, e.g. first divide indices into groups defined by partition and apply BSP on these sets
size_t hlib_clt_bytesize | ( | const hlib_clustertree_t | ct, |
int * | info | ||
) |
return size of memory in bytes used by cluster
size_t hlib_clt_depth | ( | const hlib_clustertree_t | ct, |
int * | info | ||
) |
return depth of cluster tree
void hlib_clt_free | ( | hlib_clustertree_t | ct, |
int * | info | ||
) |
free resources coupled with cluster tree ct
size_t hlib_clt_nnodes | ( | const hlib_clustertree_t | ct, |
int * | info | ||
) |
return number of nodes in cluster tree
hlib_permutation_t hlib_clt_perm_e2i | ( | hlib_clustertree_t | ct, |
int * | info | ||
) |
return mapping from external to internal (in cluster tree) ordering
hlib_permutation_t hlib_clt_perm_i2e | ( | hlib_clustertree_t | ct, |
int * | info | ||
) |
return mapping from internal (in cluster tree) to external ordering
void hlib_clt_print_ps | ( | const hlib_clustertree_t | ct, |
const char * | filename, | ||
int * | info | ||
) |
print cluster tree in PostScript format to file filename
hlib_cluster_t hlib_clt_root | ( | hlib_clustertree_t | ct, |
int * | info | ||
) |
return root node in cluster tree
hlib_real_t hlib_cmplx_abs | ( | const hlib_complex_t | a | ) |
return absolute value
hlib_complex_t hlib_cmplx_add | ( | const hlib_complex_t | a, |
const hlib_complex_t | b | ||
) |
conversion to/from C99 complex if available standard mathematical operators
hlib_complex_t hlib_cmplx_conj | ( | const hlib_complex_t | a | ) |
return conjugate value
hlib_complex_t hlib_cmplx_exp | ( | const hlib_complex_t | a | ) |
exponential function e^a
hlib_complex_t hlib_cmplx_sqrt | ( | const hlib_complex_t | a | ) |
return square root
hlib_complex_t hlib_complex | ( | const double | re, |
const double | im | ||
) |
construct hlib_complex_t variable
double* hlib_coord_bbmax | ( | const hlib_coord_t | coord, |
const size_t | i, | ||
int * | info | ||
) |
return pointer too coordinate of i'th maximal bounding box
double* hlib_coord_bbmin | ( | const hlib_coord_t | coord, |
const size_t | i, | ||
int * | info | ||
) |
return pointer too coordinate of i'th minimal bounding box
size_t hlib_coord_bytesize | ( | const hlib_coord_t | coord, |
int * | info | ||
) |
return size of memory in bytes used by coordinates
unsigned int hlib_coord_dim | ( | const hlib_coord_t | coord, |
int * | info | ||
) |
return dimension of coordinates in given coordinate set coord
void hlib_coord_free | ( | hlib_coord_t | coord, |
int * | info | ||
) |
free resources coupled with coordinates
double* hlib_coord_get | ( | const hlib_coord_t | coord, |
const size_t | i, | ||
int * | info | ||
) |
return pointer too coordinate i in given coordinate set coord
void hlib_coord_get_period | ( | hlib_coord_t | coord, |
double | period[], | ||
int * | info | ||
) |
store periodicity vector in period
int hlib_coord_has_bbox | ( | const hlib_coord_t | coord, |
int * | info | ||
) |
return 1, if coordinate set has bounding box info and 0 otherwise
hlib_coord_t hlib_coord_import | ( | const size_t | n, |
const unsigned int | dim, | ||
double ** | coord, | ||
const double * | period, | ||
int * | info | ||
) |
import given coordinates into HLIBpro
size_t hlib_coord_ncoord | ( | const hlib_coord_t | coord, |
int * | info | ||
) |
return number of coordinates in given coordinate set coord
void hlib_coord_print_vrml | ( | const hlib_coord_t | coord, |
const char * | filename, | ||
int * | info | ||
) |
print coordinates in VRML format to file filename
void hlib_coord_set_period | ( | hlib_coord_t | coord, |
const double | period[], | ||
int * | info | ||
) |
set perdiodicity of coordinates to period
double hlib_cputime | ( | ) |
return execution time of program in seconds
hlib_matrix_t hlib_hb_load_matrix | ( | const char * | filename, |
int * | info | ||
) |
read matrix from file filename
void hlib_hb_save_matrix | ( | const hlib_matrix_t | M, |
const char * | filename, | ||
int * | info | ||
) |
save matrix M to file filename
hlib_coord_t hlib_hformat_load_coord | ( | const char * | filename, |
int * | info | ||
) |
read coordinates from file filename
hlib_matrix_t hlib_hformat_load_matrix | ( | const char * | filename, |
int * | info | ||
) |
read matrix from file filename
hlib_vector_t hlib_hformat_load_vector | ( | const char * | filename, |
int * | info | ||
) |
read vector from file filename
void hlib_hformat_save_coord | ( | const hlib_coord_t | coord, |
const char * | filename, | ||
int * | info | ||
) |
save coordinates to file filename
void hlib_hformat_save_matrix | ( | const hlib_matrix_t | A, |
const char * | filename, | ||
int * | info | ||
) |
save matrix A to file filename
void hlib_hformat_save_vector | ( | const hlib_vector_t | v, |
const char * | filename, | ||
int * | info | ||
) |
save vector v to file filename
hlib_coord_t hlib_load_coord | ( | const char * | filename, |
int * | info | ||
) |
read coordinates from file filename
hlib_matrix_t hlib_load_matrix | ( | const char * | filename, |
int * | info | ||
) |
read matrix from file filename
hlib_vector_t hlib_load_vector | ( | const char * | filename, |
int * | info | ||
) |
read vector from file filename
unsigned int hlib_major_version | ( | ) |
return major version number of 𝓗𝖫𝖨𝖡𝗉𝗋𝗈
hlib_matrix_t hlib_matlab_load_matrix | ( | const char * | filename, |
const char * | matname, | ||
int * | info | ||
) |
read matrix named matname in Matlab V7 format from file filename; if matname == NULL or matname == "", the first matrix in filename will be read
hlib_vector_t hlib_matlab_load_vector | ( | const char * | filename, |
const char * | vecname, | ||
int * | info | ||
) |
read vector named vname from file filename in Matlab V7 format
void hlib_matlab_save_matrix | ( | const hlib_matrix_t | M, |
const char * | filename, | ||
const char * | matname, | ||
int * | info | ||
) |
save matrix M named matname in Matlab V7 format to file filename
void hlib_matlab_save_vector | ( | const hlib_vector_t | v, |
const char * | filename, | ||
const char * | vecname, | ||
int * | info | ||
) |
save vector named vname to file filename in Matlab V7 format
void hlib_matrix_add | ( | const hlib_real_t | alpha, |
const hlib_matrix_t | A, | ||
const hlib_real_t | beta, | ||
hlib_matrix_t | B, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
compute B := alpha A + beta B with block-wise precision acc
hlib_matrix_t hlib_matrix_assemble_block | ( | const size_t | brows, |
const size_t | bcols, | ||
hlib_matrix_t * | submat, | ||
const int | copy, | ||
int * | info | ||
) |
assemble brows × bcols block matrix out of sub blocks submat
hlib_matrix_t hlib_matrix_build_ccoeff | ( | const hlib_blockclustertree_t | bct, |
const hlib_ccoeff_t | f, | ||
void * | arg, | ||
const hlib_lrapx_t | lrapx, | ||
const hlib_acc_t | acc, | ||
const int | sym, | ||
int * | info | ||
) |
build H-matrix over block cluster tree bct of precision acc out of dense matrix given by a matrix coefficient function f by using a low-rank approximation for admissible blocks
hlib_matrix_t hlib_matrix_build_coeff | ( | const hlib_blockclustertree_t | bct, |
const hlib_coeff_t | f, | ||
void * | arg, | ||
const hlib_lrapx_t | lrapx, | ||
const hlib_acc_t | acc, | ||
const int | sym, | ||
int * | info | ||
) |
build H-matrix over block cluster tree bct of precision acc out of dense matrix given by a matrix coefficient function f by using a low-rank approximation for admissible blocks
hlib_matrix_t hlib_matrix_build_dense | ( | const hlib_blockclustertree_t | bct, |
const hlib_matrix_t | D, | ||
const hlib_lrapx_t | lrapx, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
build H-matrix over block cluster tree bct with contents defined by given dense matrix D with low-rank approximation method lrapx and block-wise accuracy acc
hlib_matrix_t hlib_matrix_build_identity | ( | const hlib_blockclustertree_t | bct, |
int * | info | ||
) |
build H-matrix over block cluster tree bct representing identity
hlib_matrix_t hlib_matrix_build_sparse | ( | const hlib_blockclustertree_t | bct, |
const hlib_matrix_t | S, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
build H-matrix over block cluster tree bct with contents defined by given sparse matrix S
size_t hlib_matrix_bytesize | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return size of matrix in bytes
void hlib_matrix_cadd | ( | const hlib_complex_t | alpha, |
const hlib_matrix_t | A, | ||
const hlib_complex_t | beta, | ||
hlib_matrix_t | B, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
compute B := alpha A + beta B with block-wise precision acc
hlib_complex_t hlib_matrix_centry_get | ( | const hlib_matrix_t | A, |
const size_t | i, | ||
const size_t | j, | ||
int * | info | ||
) |
get single entry at position (i,j) in matrix A
void hlib_matrix_cmul | ( | const hlib_complex_t | alpha, |
const hlib_matop_t | matop_A, | ||
const hlib_matrix_t | A, | ||
const hlib_matop_t | matop_B, | ||
const hlib_matrix_t | B, | ||
const hlib_complex_t | beta, | ||
hlib_matrix_t | C, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
compute C := alpha op(A) op(B) + beta C with block-wise precision acc
void hlib_matrix_cmulvec | ( | const hlib_complex_t | alpha, |
const hlib_matrix_t | A, | ||
const hlib_vector_t | x, | ||
const hlib_complex_t | beta, | ||
hlib_vector_t | y, | ||
const hlib_matop_t | matop, | ||
int * | info | ||
) |
compute y := alpha op(A) x + beta y, where op(A) is either A, A^T or A^H
hlib_vector_t hlib_matrix_col_vector | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
create vectors matching the column indexsets of matrix A
size_t hlib_matrix_cols | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return number of columns of matrix
void hlib_matrix_conjugate | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
conjugate matrix A
hlib_matrix_t hlib_matrix_copy | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return copy of matrix A
hlib_matrix_t hlib_matrix_copy_acc | ( | const hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
const int | coarsen, | ||
int * | info | ||
) |
return copy of matrix A with accuracy acc; if coarsen != 0, copy is coarsend
hlib_matrix_t hlib_matrix_copy_blockdiag | ( | const hlib_matrix_t | A, |
const unsigned int | lvl, | ||
int * | info | ||
) |
copy blockdiagonal part of first lvl levels of matrix
hlib_matrix_t hlib_matrix_copy_blockdiag_acc | ( | const hlib_matrix_t | A, |
const unsigned int | lvl, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
copy blockdiagonal part of first lvl levels of matrix with given accuracy acc
hlib_matrix_t hlib_matrix_copy_nearfield | ( | const hlib_matrix_t | A, |
const int | without_rk, | ||
int * | info | ||
) |
copy nearfield part of matrix
void hlib_matrix_copyto | ( | const hlib_matrix_t | A, |
hlib_matrix_t | B, | ||
int * | info | ||
) |
copy matrix A to B
void hlib_matrix_copyto_acc | ( | const hlib_matrix_t | A, |
hlib_matrix_t | B, | ||
const hlib_acc_t | acc, | ||
const int | coarsen, | ||
int * | info | ||
) |
copy matrix A to B with accuracy acc; if coarsen != 0, B is coarsend
void hlib_matrix_cscale | ( | const hlib_complex_t | f, |
const hlib_matrix_t | A, | ||
int * | info | ||
) |
scale matrix A by factor f
hlib_real_t hlib_matrix_entry_get | ( | const hlib_matrix_t | A, |
const size_t | i, | ||
const size_t | j, | ||
int * | info | ||
) |
get single entry at position (i,j) in matrix A
hlib_matrix_t hlib_matrix_factorise | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
factorise matrix A up to block-wise precision acc
hlib_matrix_t hlib_matrix_factorise_inv | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
same as hlib_matrix_factorise but return matrix object for evaluation of inverse of A.
void hlib_matrix_free | ( | hlib_matrix_t | A, |
int * | info | ||
) |
free resources coupled with matrix A
hlib_matrix_t hlib_matrix_import_ccrs | ( | const size_t | rows, |
const size_t | cols, | ||
const size_t | nnz, | ||
const int * | rowind, | ||
const int * | colind, | ||
const hlib_complex_t * | coeffs, | ||
const int | sym, | ||
int * | info | ||
) |
import given complex sparse matrix in CRS format into internal sparse matrix format
hlib_matrix_t hlib_matrix_import_cdense | ( | const size_t | rows, |
const size_t | cols, | ||
const hlib_complex_t * | D, | ||
const int | sym, | ||
int * | info | ||
) |
import given complex dense matrix in column major format into internal dense matrix format
hlib_matrix_t hlib_matrix_import_crs | ( | const size_t | rows, |
const size_t | cols, | ||
const size_t | nnz, | ||
const int * | rowind, | ||
const int * | colind, | ||
const hlib_real_t * | coeffs, | ||
const int | sym, | ||
int * | info | ||
) |
import given real sparse matrix in CRS format into internal sparse matrix format
hlib_matrix_t hlib_matrix_import_dense | ( | const size_t | rows, |
const size_t | cols, | ||
const hlib_real_t * | D, | ||
const int | sym, | ||
int * | info | ||
) |
import given real dense matrix in column major format into internal dense matrix format
hlib_matrix_t hlib_matrix_inv | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
compute inverse of A with block-wise precision acc; A will be overwritten with A^-1 upon exit
hlib_vector_t hlib_matrix_inv_diag | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
compute diagonal of inverse of A with local accuracy acc
int hlib_matrix_is_complex | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return 1 if matrix A is complex valued and 0 otherwise
int hlib_matrix_is_herm | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return 1 if matrix A is hermitian, e.g. A = A^H
int hlib_matrix_is_sym | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return 1 if matrix A is symmetric, e.g. A = A^T
hlib_matrix_t hlib_matrix_ldl | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
LDL factorisation of matrix A up to block-wise precision acc
hlib_matrix_t hlib_matrix_ll | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
Cholesky factorisation of matrix A up to block-wise precision acc
hlib_matrix_t hlib_matrix_lu | ( | hlib_matrix_t | A, |
const hlib_acc_t | acc, | ||
int * | info | ||
) |
LU factorise matrix A up to block-wise precision acc
void hlib_matrix_mul | ( | const hlib_real_t | alpha, |
const hlib_matop_t | matop_A, | ||
const hlib_matrix_t | A, | ||
const hlib_matop_t | matop_B, | ||
const hlib_matrix_t | B, | ||
const hlib_real_t | beta, | ||
hlib_matrix_t | C, | ||
const hlib_acc_t | acc, | ||
int * | info | ||
) |
compute C := alpha op(A) op(B) + beta C with block-wise precision acc
void hlib_matrix_mulvec | ( | const hlib_real_t | alpha, |
const hlib_matrix_t | A, | ||
const hlib_vector_t | x, | ||
const hlib_real_t | beta, | ||
hlib_vector_t | y, | ||
const hlib_matop_t | matop, | ||
int * | info | ||
) |
compute y := alpha op(A) x + beta y, where op(A) is either A or A^T
hlib_real_t hlib_matrix_norm_frobenius | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
compute Frobenius norm of matrix
hlib_real_t hlib_matrix_norm_frobenius_diff | ( | const hlib_matrix_t | A, |
const hlib_matrix_t | B, | ||
int * | info | ||
) |
compute Frobenius norm of A-B
hlib_real_t hlib_matrix_norm_inv_approx | ( | const hlib_matrix_t | A, |
const hlib_matrix_t | B, | ||
int * | info | ||
) |
compute inversion error |I-BA|_2, where B is supposed to be A^-1
hlib_real_t hlib_matrix_norm_spectral | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
compute spectral norm of matrix
hlib_real_t hlib_matrix_norm_spectral_diff | ( | const hlib_matrix_t | A, |
const hlib_matrix_t | B, | ||
int * | info | ||
) |
compute relative spectral norm of A-B, e.g. |A-B|_2/|A|_2
hlib_real_t hlib_matrix_norm_spectral_inv | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
compute spectral norm of inverse matrix
void hlib_matrix_print_ps | ( | const hlib_matrix_t | A, |
const char * | filename, | ||
const int | options, | ||
int * | info | ||
) |
print matrix A in postscript format to file filename
hlib_vector_t hlib_matrix_row_vector | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
create vectors matching the row indexsets of matrix A
size_t hlib_matrix_rows | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
return number of rows of matrix
void hlib_matrix_scale | ( | const hlib_real_t | f, |
const hlib_matrix_t | A, | ||
int * | info | ||
) |
scale matrix A by factor f
void hlib_matrix_solve_diag_left | ( | const hlib_matrix_t | D, |
const hlib_matop_t | op_D, | ||
hlib_matrix_t | A, | ||
const hlib_acc_t | acc, | ||
const int | pointwise, | ||
int * | info | ||
) |
solve op(D) B = A with diagonal D and known A
void hlib_matrix_solve_diag_right | ( | hlib_matrix_t | A, |
const hlib_matrix_t | D, | ||
const hlib_matop_t | op_D, | ||
const hlib_acc_t | acc, | ||
const int | pointwise, | ||
int * | info | ||
) |
solve B op(D) = A with diagonal D and known A
void hlib_matrix_solve_lower_left | ( | const hlib_matrix_t | L, |
hlib_matrix_t | A, | ||
const hlib_acc_t | acc, | ||
const int | pointwise, | ||
const int | unit_diag, | ||
int * | info | ||
) |
solve L B = A with lower triangular L and known A
void hlib_matrix_solve_lower_right | ( | hlib_matrix_t | A, |
const hlib_matrix_t | L, | ||
const hlib_matop_t | op_L, | ||
const hlib_acc_t | acc, | ||
const int | pointwise, | ||
const int | unit_diag, | ||
int * | info | ||
) |
solve B op(L) = A with lower triangular L and known A
void hlib_matrix_solve_upper_right | ( | hlib_matrix_t | A, |
const hlib_matrix_t | U, | ||
const hlib_acc_t | acc, | ||
const int | pointwise, | ||
const int | unit_diag, | ||
int * | info | ||
) |
solve B U = A with upper triangular U and known A
void hlib_matrix_to_complex | ( | const hlib_matrix_t | A, |
const int | force, | ||
int * | info | ||
) |
switch matrix to complex value; if force is 1, the conversion is performed independent of the current state (applies to block matrices)
void hlib_matrix_to_real | ( | const hlib_matrix_t | A, |
const int | force, | ||
int * | info | ||
) |
switch matrix to real value (if possible); if force is 1, the conversion is performed independent of the current state (applies to block matrices)
void hlib_matrix_transpose | ( | const hlib_matrix_t | A, |
int * | info | ||
) |
transpose matrix A
unsigned int hlib_minor_version | ( | ) |
return minor version number of 𝓗𝖫𝖨𝖡𝗉𝗋𝗈
hlib_matrix_t hlib_mm_load_matrix | ( | const char * | filename, |
int * | info | ||
) |
read matrix from file filename
hlib_matrix_t hlib_pltmg_load_matrix | ( | const char * | filename, |
int * | info | ||
) |
read matrix from file filename
hlib_coord_t hlib_samg_load_coord | ( | const char * | filename, |
int * | info | ||
) |
read coordinates from file filename
hlib_matrix_t hlib_samg_load_matrix | ( | const char * | filename, |
int * | info | ||
) |
read matrix A from file "\a filename" in SAMG format
hlib_vector_t hlib_samg_load_vector | ( | const char * | filename, |
int * | info | ||
) |
read vector from file filename in SAMG format
void hlib_samg_save_coord | ( | const hlib_coord_t | coord, |
const char * | filename, | ||
int * | info | ||
) |
write coordinates coord to file filename
void hlib_samg_save_matrix | ( | const hlib_matrix_t | A, |
const char * | filename, | ||
int * | info | ||
) |
save matrix A to file "\a filename" in SAMG format
void hlib_samg_save_vector | ( | const hlib_vector_t | x, |
const char * | filename, | ||
int * | info | ||
) |
save vector to file filename in SAMG format
void hlib_set_abs_eps | ( | const hlib_real_t | eps | ) |
define minimal boundary for singular values, e.g. not smaller
void hlib_set_coarsening | ( | const int | build, |
const int | arith | ||
) |
enable (1) or disable (0) coarsening during H-matrix building and arithmetic (default: build=1, arith=0)
void hlib_set_diag_scale | ( | const int | scale | ) |
enable (1) or disable (0) diagonal scaling during LU, etc. (default: off)
void hlib_set_error_fn | ( | const hlib_errorfn_t | errorfn | ) |
set call back function in case of an error
void hlib_set_n_min | ( | const unsigned int | n | ) |
set maximal leaf size in cluster trees to n
void hlib_set_nthreads | ( | const unsigned int | p | ) |
set maximal number of threads to use in H arithmetic
void hlib_set_progress_cb | ( | hlib_progressfn_t | fn, |
void * | arg | ||
) |
set callback function for progress information if fn == NULL, HLIBpro uses default progress output
void hlib_set_recompress | ( | const int | recompress | ) |
enable (1) or disable (0) recompression during H-matrix building (default: on)
void hlib_set_verbosity | ( | const unsigned int | verb | ) |
define verbosity of HLib
void hlib_set_warning_fn | ( | const hlib_errorfn_t | warnfn | ) |
set call back function in case of a warning
hlib_solver_t hlib_solver_auto | ( | int * | info | ) |
create automatic solver
hlib_solver_t hlib_solver_bicgstab | ( | int * | info | ) |
create BiCG-Stab solver
hlib_solver_t hlib_solver_cg | ( | int * | info | ) |
create CG solver
void hlib_solver_free | ( | hlib_solver_t | solver, |
int * | info | ||
) |
free solver
hlib_solver_t hlib_solver_gmres | ( | const int | restart, |
int * | info | ||
) |
create GMRES solver with restart after restart steps
hlib_solver_t hlib_solver_minres | ( | int * | info | ) |
create MINRES solver
hlib_solver_t hlib_solver_richardson | ( | int * | info | ) |
create Richardson solver
void hlib_solver_solve | ( | const hlib_solver_t | solver, |
const hlib_matrix_t | A, | ||
hlib_vector_t | x, | ||
const hlib_vector_t | b, | ||
const hlib_matrix_t | W, | ||
hlib_solve_info_t * | solve_info, | ||
int * | info | ||
) |
solve A x = b using solver with optional preconditioning, e.g. W A x = W b if W != NULL; information about the solving process is stored in solve_info
void hlib_solver_stopcrit | ( | hlib_solver_t | solver, |
const int | maxit, | ||
const hlib_real_t | abs_red, | ||
const hlib_real_t | rel_red, | ||
int * | info | ||
) |
set stopping criterion for solver
void hlib_vector_assign | ( | hlib_vector_t | y, |
const hlib_vector_t | x, | ||
int * | info | ||
) |
copy content of vector x to vector y x and y have to be of equal type
void hlib_vector_axpy | ( | const hlib_real_t | alpha, |
const hlib_vector_t | x, | ||
hlib_vector_t | y, | ||
int * | info | ||
) |
compute y := y + a x
hlib_vector_t hlib_vector_build | ( | const size_t | size, |
int * | info | ||
) |
create scalar vectors of size size initialised with 0
size_t hlib_vector_bytesize | ( | const hlib_vector_t | v, |
int * | info | ||
) |
return size of matrix in bytes
void hlib_vector_caxpy | ( | const hlib_complex_t | alpha, |
const hlib_vector_t | x, | ||
hlib_vector_t | y, | ||
int * | info | ||
) |
compute y := y + a x
hlib_vector_t hlib_vector_cbuild | ( | const size_t | size, |
int * | info | ||
) |
create scalar vectors of size size initialised with 0
hlib_complex_t hlib_vector_centry_get | ( | const hlib_vector_t | x, |
const size_t | i, | ||
int * | info | ||
) |
get single entry at position i in vector x
void hlib_vector_centry_set | ( | const hlib_vector_t | x, |
const size_t | i, | ||
const hlib_complex_t | f, | ||
int * | info | ||
) |
set single entry at position i in vector x
void hlib_vector_cexport | ( | const hlib_vector_t | v, |
hlib_complex_t * | arr, | ||
int * | info | ||
) |
copy data from vector v into given C arrays
void hlib_vector_cfill | ( | hlib_vector_t | x, |
const hlib_complex_t | f, | ||
int * | info | ||
) |
fill vector x with constant value f
void hlib_vector_cimport | ( | hlib_vector_t | v, |
const hlib_complex_t * | arr, | ||
int * | info | ||
) |
copy data from C array arr into vector v
hlib_vector_t hlib_vector_copy | ( | const hlib_vector_t | v, |
int * | info | ||
) |
copy vectormatrix
void hlib_vector_cscale | ( | hlib_vector_t | x, |
const hlib_complex_t | f, | ||
int * | info | ||
) |
scale vector x by f
hlib_complex_t hlib_vector_dot | ( | const hlib_vector_t | x, |
const hlib_vector_t | y, | ||
int * | info | ||
) |
compute dot product <x,y> = x^H * y
hlib_real_t hlib_vector_entry_get | ( | const hlib_vector_t | x, |
const size_t | i, | ||
int * | info | ||
) |
get single entry at position i in vector x
void hlib_vector_entry_set | ( | const hlib_vector_t | x, |
const size_t | i, | ||
const hlib_real_t | f, | ||
int * | info | ||
) |
set single entry at position i in vector x
void hlib_vector_export | ( | const hlib_vector_t | v, |
hlib_real_t * | arr, | ||
int * | info | ||
) |
copy data from vector v into given C arrays
void hlib_vector_fft | ( | const hlib_vector_t | v, |
int * | info | ||
) |
compute FFT of vector v inplace
void hlib_vector_fill | ( | hlib_vector_t | x, |
const hlib_real_t | f, | ||
int * | info | ||
) |
fill vector x with constant value f
void hlib_vector_fill_rand | ( | hlib_vector_t | x, |
int * | info | ||
) |
fill vector x with random values
void hlib_vector_free | ( | hlib_vector_t | v, |
int * | info | ||
) |
free resources coupled with vector v
void hlib_vector_ifft | ( | const hlib_vector_t | v, |
int * | info | ||
) |
compute inverse FFT of vector v inplace
void hlib_vector_import | ( | hlib_vector_t | v, |
const hlib_real_t * | arr, | ||
int * | info | ||
) |
copy data from C array arr into vector v
hlib_real_t hlib_vector_norm2 | ( | const hlib_vector_t | x, |
int * | info | ||
) |
compute euklidean norm of x
hlib_real_t hlib_vector_norm_inf | ( | const hlib_vector_t | x, |
int * | info | ||
) |
compute infinity norm of x
hlib_vector_t hlib_vector_restrict_im | ( | const hlib_vector_t | v, |
int * | info | ||
) |
return vector containing only imaginary parts * of the entries of the given vector
hlib_vector_t hlib_vector_restrict_re | ( | const hlib_vector_t | v, |
int * | info | ||
) |
return vector containing only real parts * of the entries of the given vector
void hlib_vector_scale | ( | hlib_vector_t | x, |
const hlib_real_t | f, | ||
int * | info | ||
) |
scale vector x by f
size_t hlib_vector_size | ( | const hlib_vector_t | v, |
int * | info | ||
) |
return size of vectors
double hlib_walltime | ( | ) |
return current walltime in seconds