:orphan: # MatCreateH2OpusFromKernel Creates a `MATH2OPUS` from a user-supplied kernel. ## Synopsis ``` PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, MatH2OpusKernel kernel, void *kernelctx, PetscReal eta, PetscInt leafsize, PetscInt basisord, Mat *nA) ``` ## Input Parameters - ***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 -*** coordinates of the points - ***cdist -*** whether or not coordinates are distributed - ***kernel -*** computational kernel (or `NULL`) - ***kernelctx -*** kernel context - ***eta -*** admissibility condition tolerance - ***leafsize -*** leaf size in cluster tree - ***basisord -*** approximation order for Chebychev interpolation of low-rank blocks ## Output Parameter - ***nA -*** matrix ## Options Database Keys - ***-mat_h2opus_leafsize <`PetscInt`> -*** Leaf size of cluster tree - ***-mat_h2opus_eta <`PetscReal`> -*** Admissibility condition tolerance - ***-mat_h2opus_order <`PetscInt`> -*** Chebychev approximation order - ***-mat_h2opus_normsamples <`PetscInt`> -*** Maximum number of samples to be used when estimating norms ## See Also [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()` ## Level intermediate ## Location src/mat/impls/h2opus/cuda/math2opus.cu ## Examples src/ksp/ksp/tutorials/ex21.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/h2opus/cuda/math2opus.cu) [Index of all Mat routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)