#include "petscmat.h" PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt spacedim,const PetscReal coords_target[],const PetscReal coords_source[],MatHtoolKernel kernel,void *kernelctx,Mat *B)
comm | - MPI communicator | |
m | - number of local rows (or PETSC_DECIDE to have calculated if M is given) | |
n | - number of local columns (or PETSC_DECIDE to have calculated if N is given) | |
M | - number of global rows (or PETSC_DETERMINE to have calculated if m is given) | |
N | - number of global columns (or PETSC_DETERMINE to have calculated if n is given) | |
spacedim | - dimension of the space coordinates | |
coords_target | - coordinates of the target | |
coords_source | - coordinates of the source | |
kernel | - computational kernel (or NULL) | |
kernelctx | - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) |
B | - matrix |
-mat_htool_min_cluster_size <PetscInt> | - minimal leaf size in cluster tree | |
-mat_htool_max_block_size <PetscInt> | - maximal number of coefficients in a dense block | |
-mat_htool_epsilon <PetscReal> | - relative error in Frobenius norm when approximating a block | |
-mat_htool_eta <PetscReal> | - admissibility condition tolerance | |
-mat_htool_min_target_depth <PetscInt> | - minimal cluster tree depth associated with the rows | |
-mat_htool_min_source_depth <PetscInt> | - minimal cluster tree depth associated with the columns | |
-mat_htool_compressor <sympartialACA, fullACA, SVD> | - type of compression | |
-mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> | - type of clustering |